ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2251
Committed: Sun May 29 21:16:25 2005 UTC (19 years, 1 month ago) by chrisfen
File size: 30578 byte(s)
Log Message:
Removed balance from the Darkside (files)

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)
311
312 logical, intent(in) :: do_pot
313
314 integer, intent(in) :: atom1, atom2
315 integer :: localError
316
317 real(kind=dp), intent(in) :: rij, r2, sw
318 real(kind=dp), intent(in), dimension(3) :: d
319 real(kind=dp), intent(inout) :: vpair
320 real(kind=dp), intent(inout), dimension(3) :: fpair
321
322 real( kind = dp ) :: pot
323 real( kind = dp ), dimension(9,nLocal) :: eFrame
324 real( kind = dp ), dimension(3,nLocal) :: f
325 real( kind = dp ), dimension(3,nLocal) :: t
326
327 real (kind = dp), dimension(3) :: ux_i, uy_i, uz_i
328 real (kind = dp), dimension(3) :: ux_j, uy_j, uz_j
329 real (kind = dp), dimension(3) :: dudux_i, duduy_i, duduz_i
330 real (kind = dp), dimension(3) :: dudux_j, duduy_j, duduz_j
331
332 logical :: i_is_Charge, i_is_Dipole, i_is_SplitDipole, i_is_Quadrupole
333 logical :: j_is_Charge, j_is_Dipole, j_is_SplitDipole, j_is_Quadrupole
334 logical :: i_is_Tap, j_is_Tap
335 integer :: me1, me2, id1, id2
336 real (kind=dp) :: q_i, q_j, mu_i, mu_j, d_i, d_j
337 real (kind=dp) :: qxx_i, qyy_i, qzz_i
338 real (kind=dp) :: qxx_j, qyy_j, qzz_j
339 real (kind=dp) :: cx_i, cy_i, cz_i
340 real (kind=dp) :: cx_j, cy_j, cz_j
341 real (kind=dp) :: cx2, cy2, cz2
342 real (kind=dp) :: ct_i, ct_j, ct_ij, a1
343 real (kind=dp) :: riji, ri, ri2, ri3, ri4
344 real (kind=dp) :: pref, vterm, epot, dudr
345 real (kind=dp) :: xhat, yhat, zhat
346 real (kind=dp) :: dudx, dudy, dudz
347 real (kind=dp) :: scale, sc2, bigR, switcher, dswitcher
348
349 if (.not.allocated(ElectrostaticMap)) then
350 call handleError("electrostatic", "no ElectrostaticMap was present before first call of do_electrostatic_pair!")
351 return
352 end if
353
354 #ifdef IS_MPI
355 me1 = atid_Row(atom1)
356 me2 = atid_Col(atom2)
357 #else
358 me1 = atid(atom1)
359 me2 = atid(atom2)
360 #endif
361
362 !! some variables we'll need independent of electrostatic type:
363
364 riji = 1.0d0 / rij
365
366 xhat = d(1) * riji
367 yhat = d(2) * riji
368 zhat = d(3) * riji
369
370 !! logicals
371 i_is_Charge = ElectrostaticMap(me1)%is_Charge
372 i_is_Dipole = ElectrostaticMap(me1)%is_Dipole
373 i_is_SplitDipole = ElectrostaticMap(me1)%is_SplitDipole
374 i_is_Quadrupole = ElectrostaticMap(me1)%is_Quadrupole
375 i_is_Tap = ElectrostaticMap(me1)%is_Tap
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 j_is_Tap = ElectrostaticMap(me2)%is_Tap
382
383 if (i_is_Charge) then
384 q_i = ElectrostaticMap(me1)%charge
385 endif
386
387 if (i_is_Dipole) then
388 mu_i = ElectrostaticMap(me1)%dipole_moment
389 #ifdef IS_MPI
390 uz_i(1) = eFrame_Row(3,atom1)
391 uz_i(2) = eFrame_Row(6,atom1)
392 uz_i(3) = eFrame_Row(9,atom1)
393 #else
394 uz_i(1) = eFrame(3,atom1)
395 uz_i(2) = eFrame(6,atom1)
396 uz_i(3) = eFrame(9,atom1)
397 #endif
398 ct_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
399
400 if (i_is_SplitDipole) then
401 d_i = ElectrostaticMap(me1)%split_dipole_distance
402 endif
403
404 endif
405
406 if (i_is_Quadrupole) then
407 qxx_i = ElectrostaticMap(me1)%quadrupole_moments(1)
408 qyy_i = ElectrostaticMap(me1)%quadrupole_moments(2)
409 qzz_i = ElectrostaticMap(me1)%quadrupole_moments(3)
410 #ifdef IS_MPI
411 ux_i(1) = eFrame_Row(1,atom1)
412 ux_i(2) = eFrame_Row(4,atom1)
413 ux_i(3) = eFrame_Row(7,atom1)
414 uy_i(1) = eFrame_Row(2,atom1)
415 uy_i(2) = eFrame_Row(5,atom1)
416 uy_i(3) = eFrame_Row(8,atom1)
417 uz_i(1) = eFrame_Row(3,atom1)
418 uz_i(2) = eFrame_Row(6,atom1)
419 uz_i(3) = eFrame_Row(9,atom1)
420 #else
421 ux_i(1) = eFrame(1,atom1)
422 ux_i(2) = eFrame(4,atom1)
423 ux_i(3) = eFrame(7,atom1)
424 uy_i(1) = eFrame(2,atom1)
425 uy_i(2) = eFrame(5,atom1)
426 uy_i(3) = eFrame(8,atom1)
427 uz_i(1) = eFrame(3,atom1)
428 uz_i(2) = eFrame(6,atom1)
429 uz_i(3) = eFrame(9,atom1)
430 #endif
431 cx_i = ux_i(1)*xhat + ux_i(2)*yhat + ux_i(3)*zhat
432 cy_i = uy_i(1)*xhat + uy_i(2)*yhat + uy_i(3)*zhat
433 cz_i = uz_i(1)*xhat + uz_i(2)*yhat + uz_i(3)*zhat
434 endif
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)*xhat + uz_j(2)*yhat + uz_j(3)*zhat
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 !!$ switcher = 1.0d0
489 !!$ dswitcher = 0.0d0
490 !!$ ebalance = 0.0d0
491 !!$ ! weaken the dipole interaction at close range for TAP water
492 !!$ if (j_is_Tap .and. i_is_Tap) then
493 !!$ call calc_switch(rij, mu_i, switcher, dswitcher)
494 !!$ endif
495
496 epot = 0.0_dp
497 dudx = 0.0_dp
498 dudy = 0.0_dp
499 dudz = 0.0_dp
500
501 dudux_i = 0.0_dp
502 duduy_i = 0.0_dp
503 duduz_i = 0.0_dp
504
505 dudux_j = 0.0_dp
506 duduy_j = 0.0_dp
507 duduz_j = 0.0_dp
508
509 if (i_is_Charge) then
510
511 if (j_is_Charge) then
512
513 vterm = pre11 * q_i * q_j * riji
514 vpair = vpair + vterm
515 epot = epot + sw*vterm
516
517 dudr = - sw * vterm * riji
518
519 dudx = dudx + dudr * xhat
520 dudy = dudy + dudr * yhat
521 dudz = dudz + dudr * zhat
522
523 endif
524
525 if (j_is_Dipole) then
526
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
536 ri2 = ri * ri
537 ri3 = ri2 * ri
538 sc2 = scale * scale
539
540 pref = pre12 * q_i * mu_j
541 vterm = - pref * ct_j * ri2 * scale
542 vpair = vpair + vterm
543 epot = epot + sw * vterm
544
545 !! this has a + sign in the () because the rij vector is
546 !! r_j - r_i and the charge-dipole potential takes the origin
547 !! as the point dipole, which is atom j in this case.
548
549 dudx = dudx - pref * sw * ri3 * ( uz_j(1) - 3.0d0*ct_j*xhat*sc2)
550 dudy = dudy - pref * sw * ri3 * ( uz_j(2) - 3.0d0*ct_j*yhat*sc2)
551 dudz = dudz - pref * sw * ri3 * ( uz_j(3) - 3.0d0*ct_j*zhat*sc2)
552
553 duduz_j(1) = duduz_j(1) - pref * sw * ri2 * xhat * scale
554 duduz_j(2) = duduz_j(2) - pref * sw * ri2 * yhat * scale
555 duduz_j(3) = duduz_j(3) - pref * sw * ri2 * zhat * scale
556
557 endif
558
559 if (j_is_Quadrupole) then
560 ri2 = riji * riji
561 ri3 = ri2 * riji
562 ri4 = ri2 * ri2
563 cx2 = cx_j * cx_j
564 cy2 = cy_j * cy_j
565 cz2 = cz_j * cz_j
566
567
568 pref = pre14 * q_i / 3.0_dp
569 vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
570 qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
571 qzz_j * (3.0_dp*cz2 - 1.0_dp))
572 vpair = vpair + vterm
573 epot = epot + sw * vterm
574
575 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + pref * sw * ri4 * ( &
576 qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + &
577 qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + &
578 qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) )
579 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + pref * sw * ri4 * ( &
580 qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + &
581 qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + &
582 qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) )
583 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + pref * sw * ri4 * ( &
584 qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + &
585 qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + &
586 qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) )
587
588 dudux_j(1) = dudux_j(1) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*xhat)
589 dudux_j(2) = dudux_j(2) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*yhat)
590 dudux_j(3) = dudux_j(3) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*zhat)
591
592 duduy_j(1) = duduy_j(1) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*xhat)
593 duduy_j(2) = duduy_j(2) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*yhat)
594 duduy_j(3) = duduy_j(3) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*zhat)
595
596 duduz_j(1) = duduz_j(1) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*xhat)
597 duduz_j(2) = duduz_j(2) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*yhat)
598 duduz_j(3) = duduz_j(3) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*zhat)
599 endif
600
601 endif
602
603 if (i_is_Dipole) then
604
605 if (j_is_Charge) then
606
607 if (i_is_SplitDipole) then
608 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
609 ri = 1.0_dp / BigR
610 scale = rij * ri
611 else
612 ri = riji
613 scale = 1.0_dp
614 endif
615
616 ri2 = ri * ri
617 ri3 = ri2 * ri
618 sc2 = scale * scale
619
620 pref = pre12 * q_j * mu_i
621 vterm = pref * ct_i * ri2 * scale
622 vpair = vpair + vterm
623 epot = epot + sw * vterm
624
625 dudx = dudx + pref * sw * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2)
626 dudy = dudy + pref * sw * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2)
627 dudz = dudz + pref * sw * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2)
628
629 duduz_i(1) = duduz_i(1) + pref * sw * ri2 * xhat * scale
630 duduz_i(2) = duduz_i(2) + pref * sw * ri2 * yhat * scale
631 duduz_i(3) = duduz_i(3) + pref * sw * ri2 * zhat * scale
632 endif
633
634 if (j_is_Dipole) then
635
636 if (i_is_SplitDipole) then
637 if (j_is_SplitDipole) then
638 BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
639 else
640 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
641 endif
642 ri = 1.0_dp / BigR
643 scale = rij * ri
644 else
645 if (j_is_SplitDipole) then
646 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
647 ri = 1.0_dp / BigR
648 scale = rij * ri
649 else
650 ri = riji
651 scale = 1.0_dp
652 endif
653 endif
654
655 ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
656
657 ri2 = ri * ri
658 ri3 = ri2 * ri
659 ri4 = ri2 * ri2
660 sc2 = scale * scale
661
662 pref = pre22 * mu_i * mu_j
663 vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
664 vpair = vpair + vterm
665 epot = epot + sw * vterm
666
667 a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
668
669 dudx = dudx + pref*sw*3.0d0*ri4*scale &
670 *(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1))
671 dudy = dudy + pref*sw*3.0d0*ri4*scale &
672 *(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2))
673 dudz = dudz + pref*sw*3.0d0*ri4*scale &
674 *(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3))
675
676 duduz_i(1) = duduz_i(1) + pref*sw*ri3 &
677 *(uz_j(1) - 3.0d0*ct_j*xhat*sc2)
678 duduz_i(2) = duduz_i(2) + pref*sw*ri3 &
679 *(uz_j(2) - 3.0d0*ct_j*yhat*sc2)
680 duduz_i(3) = duduz_i(3) + pref*sw*ri3 &
681 *(uz_j(3) - 3.0d0*ct_j*zhat*sc2)
682
683 duduz_j(1) = duduz_j(1) + pref*sw*ri3 &
684 *(uz_i(1) - 3.0d0*ct_i*xhat*sc2)
685 duduz_j(2) = duduz_j(2) + pref*sw*ri3 &
686 *(uz_i(2) - 3.0d0*ct_i*yhat*sc2)
687 duduz_j(3) = duduz_j(3) + pref*sw*ri3 &
688 *(uz_i(3) - 3.0d0*ct_i*zhat*sc2)
689 endif
690
691 endif
692
693 if (i_is_Quadrupole) then
694 if (j_is_Charge) then
695
696 ri2 = riji * riji
697 ri3 = ri2 * riji
698 ri4 = ri2 * ri2
699 cx2 = cx_i * cx_i
700 cy2 = cy_i * cy_i
701 cz2 = cz_i * cz_i
702
703 pref = pre14 * q_j / 3.0_dp
704 vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
705 qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
706 qzz_i * (3.0_dp*cz2 - 1.0_dp))
707 vpair = vpair + vterm
708 epot = epot + sw * vterm
709
710 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + pref * sw * ri4 * ( &
711 qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + &
712 qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + &
713 qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) )
714 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + pref * sw * ri4 * ( &
715 qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + &
716 qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + &
717 qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) )
718 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + pref * sw * ri4 * ( &
719 qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + &
720 qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + &
721 qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) )
722
723 dudux_i(1) = dudux_i(1) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*xhat)
724 dudux_i(2) = dudux_i(2) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*yhat)
725 dudux_i(3) = dudux_i(3) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*zhat)
726
727 duduy_i(1) = duduy_i(1) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*xhat)
728 duduy_i(2) = duduy_i(2) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*yhat)
729 duduy_i(3) = duduy_i(3) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*zhat)
730
731 duduz_i(1) = duduz_i(1) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*xhat)
732 duduz_i(2) = duduz_i(2) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*yhat)
733 duduz_i(3) = duduz_i(3) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*zhat)
734 endif
735 endif
736
737
738 if (do_pot) then
739 #ifdef IS_MPI
740 pot_row(atom1) = pot_row(atom1) + 0.5d0*epot
741 pot_col(atom2) = pot_col(atom2) + 0.5d0*epot
742 #else
743 pot = pot + epot
744 #endif
745 endif
746
747 #ifdef IS_MPI
748 f_Row(1,atom1) = f_Row(1,atom1) + dudx
749 f_Row(2,atom1) = f_Row(2,atom1) + dudy
750 f_Row(3,atom1) = f_Row(3,atom1) + dudz
751
752 f_Col(1,atom2) = f_Col(1,atom2) - dudx
753 f_Col(2,atom2) = f_Col(2,atom2) - dudy
754 f_Col(3,atom2) = f_Col(3,atom2) - dudz
755
756 if (i_is_Dipole .or. i_is_Quadrupole) then
757 t_Row(1,atom1)=t_Row(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
758 t_Row(2,atom1)=t_Row(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
759 t_Row(3,atom1)=t_Row(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
760 endif
761 if (i_is_Quadrupole) then
762 t_Row(1,atom1)=t_Row(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
763 t_Row(2,atom1)=t_Row(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
764 t_Row(3,atom1)=t_Row(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
765
766 t_Row(1,atom1)=t_Row(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
767 t_Row(2,atom1)=t_Row(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
768 t_Row(3,atom1)=t_Row(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
769 endif
770
771 if (j_is_Dipole .or. j_is_Quadrupole) then
772 t_Col(1,atom2)=t_Col(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
773 t_Col(2,atom2)=t_Col(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
774 t_Col(3,atom2)=t_Col(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
775 endif
776 if (j_is_Quadrupole) then
777 t_Col(1,atom2)=t_Col(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
778 t_Col(2,atom2)=t_Col(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
779 t_Col(3,atom2)=t_Col(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
780
781 t_Col(1,atom2)=t_Col(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
782 t_Col(2,atom2)=t_Col(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
783 t_Col(3,atom2)=t_Col(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
784 endif
785
786 #else
787 f(1,atom1) = f(1,atom1) + dudx
788 f(2,atom1) = f(2,atom1) + dudy
789 f(3,atom1) = f(3,atom1) + dudz
790
791 f(1,atom2) = f(1,atom2) - dudx
792 f(2,atom2) = f(2,atom2) - dudy
793 f(3,atom2) = f(3,atom2) - dudz
794
795 if (i_is_Dipole .or. i_is_Quadrupole) then
796 t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
797 t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
798 t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
799 endif
800 if (i_is_Quadrupole) then
801 t(1,atom1)=t(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
802 t(2,atom1)=t(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
803 t(3,atom1)=t(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
804
805 t(1,atom1)=t(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
806 t(2,atom1)=t(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
807 t(3,atom1)=t(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
808 endif
809
810 if (j_is_Dipole .or. j_is_Quadrupole) then
811 t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
812 t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
813 t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
814 endif
815 if (j_is_Quadrupole) then
816 t(1,atom2)=t(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
817 t(2,atom2)=t(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
818 t(3,atom2)=t(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
819
820 t(1,atom2)=t(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
821 t(2,atom2)=t(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
822 t(3,atom2)=t(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
823 endif
824
825 #endif
826
827 #ifdef IS_MPI
828 id1 = AtomRowToGlobal(atom1)
829 id2 = AtomColToGlobal(atom2)
830 #else
831 id1 = atom1
832 id2 = atom2
833 #endif
834
835 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
836
837 fpair(1) = fpair(1) + dudx
838 fpair(2) = fpair(2) + dudy
839 fpair(3) = fpair(3) + dudz
840
841 endif
842
843 return
844 end subroutine doElectrostaticPair
845
846 !! calculates the switching functions and their derivatives for a given
847 subroutine calc_switch(r, mu, scale, dscale)
848
849 real (kind=dp), intent(in) :: r, mu
850 real (kind=dp), intent(inout) :: scale, dscale
851 real (kind=dp) :: rl, ru, mulow, minRatio, temp, scaleVal
852
853 ! distances must be in angstroms
854 rl = 2.75d0
855 ru = 3.75d0
856 mulow = 0.0d0 !3.3856d0 ! 1.84 * 1.84
857 minRatio = mulow / (mu*mu)
858 scaleVal = 1.0d0 - minRatio
859
860 if (r.lt.rl) then
861 scale = minRatio
862 dscale = 0.0d0
863 elseif (r.gt.ru) then
864 scale = 1.0d0
865 dscale = 0.0d0
866 else
867 scale = 1.0d0 - scaleVal*((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2) &
868 / ((ru - rl)**3)
869 dscale = -scaleVal * 6.0d0 * (r-ru)*(r-rl)/((ru - rl)**3)
870 endif
871
872 return
873 end subroutine calc_switch
874
875 subroutine destroyElectrostaticTypes()
876
877 if(allocated(ElectrostaticMap)) deallocate(ElectrostaticMap)
878
879 end subroutine destroyElectrostaticTypes
880
881 end module electrostatic_module