ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/electrostatic.F90
Revision: 2229
Committed: Tue May 17 22:35:01 2005 UTC (19 years, 2 months ago) by chrisfen
File size: 30603 byte(s)
Log Message:
Modifications to tap.  Also correcting changes to the previous merge that were not caught

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 ! if (j_is_Tap .and. i_is_Tap) then
491 ! call calc_switch(rij, mu_i, switcher, dswitcher)
492 ! endif
493
494 epot = 0.0_dp
495 dudx = 0.0_dp
496 dudy = 0.0_dp
497 dudz = 0.0_dp
498
499 dudux_i = 0.0_dp
500 duduy_i = 0.0_dp
501 duduz_i = 0.0_dp
502
503 dudux_j = 0.0_dp
504 duduy_j = 0.0_dp
505 duduz_j = 0.0_dp
506
507 if (i_is_Charge) then
508
509 if (j_is_Charge) then
510
511 vterm = pre11 * q_i * q_j * riji
512 vpair = vpair + vterm
513 epot = epot + sw*vterm
514
515 dudr = - sw * vterm * riji
516
517 dudx = dudx + dudr * xhat
518 dudy = dudy + dudr * yhat
519 dudz = dudz + dudr * zhat
520
521 endif
522
523 if (j_is_Dipole) then
524
525 if (j_is_SplitDipole) then
526 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
527 ri = 1.0_dp / BigR
528 scale = rij * ri
529 else
530 ri = riji
531 scale = 1.0_dp
532 endif
533
534 ri2 = ri * ri
535 ri3 = ri2 * ri
536 sc2 = scale * scale
537
538 pref = pre12 * q_i * mu_j
539 vterm = - pref * ct_j * ri2 * scale
540 vpair = vpair + vterm
541 epot = epot + sw * vterm
542
543 !! this has a + sign in the () because the rij vector is
544 !! r_j - r_i and the charge-dipole potential takes the origin
545 !! as the point dipole, which is atom j in this case.
546
547 dudx = dudx - pref * sw * ri3 * ( uz_j(1) - 3.0d0*ct_j*xhat*sc2)
548 dudy = dudy - pref * sw * ri3 * ( uz_j(2) - 3.0d0*ct_j*yhat*sc2)
549 dudz = dudz - pref * sw * ri3 * ( uz_j(3) - 3.0d0*ct_j*zhat*sc2)
550
551 duduz_j(1) = duduz_j(1) - pref * sw * ri2 * xhat * scale
552 duduz_j(2) = duduz_j(2) - pref * sw * ri2 * yhat * scale
553 duduz_j(3) = duduz_j(3) - pref * sw * ri2 * zhat * scale
554
555 endif
556
557 if (j_is_Quadrupole) then
558 ri2 = riji * riji
559 ri3 = ri2 * riji
560 ri4 = ri2 * ri2
561 cx2 = cx_j * cx_j
562 cy2 = cy_j * cy_j
563 cz2 = cz_j * cz_j
564
565
566 pref = pre14 * q_i / 3.0_dp
567 vterm = pref * ri3 * (qxx_j * (3.0_dp*cx2 - 1.0_dp) + &
568 qyy_j * (3.0_dp*cy2 - 1.0_dp) + &
569 qzz_j * (3.0_dp*cz2 - 1.0_dp))
570 vpair = vpair + vterm
571 epot = epot + sw * vterm
572
573 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + pref * sw * ri4 * ( &
574 qxx_j*(6.0_dp*cx_j*ux_j(1) - 2.0_dp*xhat) + &
575 qyy_j*(6.0_dp*cy_j*uy_j(1) - 2.0_dp*xhat) + &
576 qzz_j*(6.0_dp*cz_j*uz_j(1) - 2.0_dp*xhat) )
577 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + pref * sw * ri4 * ( &
578 qxx_j*(6.0_dp*cx_j*ux_j(2) - 2.0_dp*yhat) + &
579 qyy_j*(6.0_dp*cy_j*uy_j(2) - 2.0_dp*yhat) + &
580 qzz_j*(6.0_dp*cz_j*uz_j(2) - 2.0_dp*yhat) )
581 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + pref * sw * ri4 * ( &
582 qxx_j*(6.0_dp*cx_j*ux_j(3) - 2.0_dp*zhat) + &
583 qyy_j*(6.0_dp*cy_j*uy_j(3) - 2.0_dp*zhat) + &
584 qzz_j*(6.0_dp*cz_j*uz_j(3) - 2.0_dp*zhat) )
585
586 dudux_j(1) = dudux_j(1) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*xhat)
587 dudux_j(2) = dudux_j(2) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*yhat)
588 dudux_j(3) = dudux_j(3) + pref * sw * ri3 * (qxx_j*6.0_dp*cx_j*zhat)
589
590 duduy_j(1) = duduy_j(1) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*xhat)
591 duduy_j(2) = duduy_j(2) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*yhat)
592 duduy_j(3) = duduy_j(3) + pref * sw * ri3 * (qyy_j*6.0_dp*cy_j*zhat)
593
594 duduz_j(1) = duduz_j(1) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*xhat)
595 duduz_j(2) = duduz_j(2) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*yhat)
596 duduz_j(3) = duduz_j(3) + pref * sw * ri3 * (qzz_j*6.0_dp*cz_j*zhat)
597 endif
598
599 endif
600
601 if (i_is_Dipole) then
602
603 if (j_is_Charge) then
604
605 if (i_is_SplitDipole) then
606 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
607 ri = 1.0_dp / BigR
608 scale = rij * ri
609 else
610 ri = riji
611 scale = 1.0_dp
612 endif
613
614 ri2 = ri * ri
615 ri3 = ri2 * ri
616 sc2 = scale * scale
617
618 pref = pre12 * q_j * mu_i
619 vterm = pref * ct_i * ri2 * scale
620 vpair = vpair + vterm
621 epot = epot + sw * vterm
622
623 dudx = dudx + pref * sw * ri3 * ( uz_i(1) - 3.0d0 * ct_i * xhat*sc2)
624 dudy = dudy + pref * sw * ri3 * ( uz_i(2) - 3.0d0 * ct_i * yhat*sc2)
625 dudz = dudz + pref * sw * ri3 * ( uz_i(3) - 3.0d0 * ct_i * zhat*sc2)
626
627 duduz_i(1) = duduz_i(1) + pref * sw * ri2 * xhat * scale
628 duduz_i(2) = duduz_i(2) + pref * sw * ri2 * yhat * scale
629 duduz_i(3) = duduz_i(3) + pref * sw * ri2 * zhat * scale
630 endif
631
632 if (j_is_Dipole) then
633
634 if (i_is_SplitDipole) then
635 if (j_is_SplitDipole) then
636 BigR = sqrt(r2 + 0.25_dp * d_i * d_i + 0.25_dp * d_j * d_j)
637 else
638 BigR = sqrt(r2 + 0.25_dp * d_i * d_i)
639 endif
640 ri = 1.0_dp / BigR
641 scale = rij * ri
642 else
643 if (j_is_SplitDipole) then
644 BigR = sqrt(r2 + 0.25_dp * d_j * d_j)
645 ri = 1.0_dp / BigR
646 scale = rij * ri
647 else
648 ri = riji
649 scale = 1.0_dp
650 endif
651 endif
652
653 ct_ij = uz_i(1)*uz_j(1) + uz_i(2)*uz_j(2) + uz_i(3)*uz_j(3)
654
655 ri2 = ri * ri
656 ri3 = ri2 * ri
657 ri4 = ri2 * ri2
658 sc2 = scale * scale
659
660 pref = pre22 * mu_i * mu_j
661 vterm = pref * ri3 * (ct_ij - 3.0d0 * ct_i * ct_j * sc2)
662 vpair = vpair + vterm
663 epot = epot + sw * vterm
664
665 a1 = 5.0d0 * ct_i * ct_j * sc2 - ct_ij
666
667 dudx = dudx + switcher*pref*sw*3.0d0*ri4*scale &
668 *(a1*xhat-ct_i*uz_j(1)-ct_j*uz_i(1)) + vterm*dswitcher
669 dudy = dudy + switcher*pref*sw*3.0d0*ri4*scale &
670 *(a1*yhat-ct_i*uz_j(2)-ct_j*uz_i(2)) + vterm*dswitcher
671 dudz = dudz + switcher*pref*sw*3.0d0*ri4*scale &
672 *(a1*zhat-ct_i*uz_j(3)-ct_j*uz_i(3)) + vterm*dswitcher
673
674 duduz_i(1) = duduz_i(1) + switcher*pref*sw*ri3 &
675 *(uz_j(1) - 3.0d0*ct_j*xhat*sc2)
676 duduz_i(2) = duduz_i(2) + switcher*pref*sw*ri3 &
677 *(uz_j(2) - 3.0d0*ct_j*yhat*sc2)
678 duduz_i(3) = duduz_i(3) + switcher*pref*sw*ri3 &
679 *(uz_j(3) - 3.0d0*ct_j*zhat*sc2)
680
681 duduz_j(1) = duduz_j(1) + switcher*pref*sw*ri3 &
682 *(uz_i(1) - 3.0d0*ct_i*xhat*sc2)
683 duduz_j(2) = duduz_j(2) + switcher*pref*sw*ri3 &
684 *(uz_i(2) - 3.0d0*ct_i*yhat*sc2)
685 duduz_j(3) = duduz_j(3) + switcher*pref*sw*ri3 &
686 *(uz_i(3) - 3.0d0*ct_i*zhat*sc2)
687 endif
688
689 endif
690
691 if (i_is_Quadrupole) then
692 if (j_is_Charge) then
693
694 ri2 = riji * riji
695 ri3 = ri2 * riji
696 ri4 = ri2 * ri2
697 cx2 = cx_i * cx_i
698 cy2 = cy_i * cy_i
699 cz2 = cz_i * cz_i
700
701 pref = pre14 * q_j / 3.0_dp
702 vterm = pref * ri3 * (qxx_i * (3.0_dp*cx2 - 1.0_dp) + &
703 qyy_i * (3.0_dp*cy2 - 1.0_dp) + &
704 qzz_i * (3.0_dp*cz2 - 1.0_dp))
705 vpair = vpair + vterm
706 epot = epot + sw * vterm
707
708 dudx = dudx - 5.0_dp*sw*vterm*riji*xhat + pref * sw * ri4 * ( &
709 qxx_i*(6.0_dp*cx_i*ux_i(1) - 2.0_dp*xhat) + &
710 qyy_i*(6.0_dp*cy_i*uy_i(1) - 2.0_dp*xhat) + &
711 qzz_i*(6.0_dp*cz_i*uz_i(1) - 2.0_dp*xhat) )
712 dudy = dudy - 5.0_dp*sw*vterm*riji*yhat + pref * sw * ri4 * ( &
713 qxx_i*(6.0_dp*cx_i*ux_i(2) - 2.0_dp*yhat) + &
714 qyy_i*(6.0_dp*cy_i*uy_i(2) - 2.0_dp*yhat) + &
715 qzz_i*(6.0_dp*cz_i*uz_i(2) - 2.0_dp*yhat) )
716 dudz = dudz - 5.0_dp*sw*vterm*riji*zhat + pref * sw * ri4 * ( &
717 qxx_i*(6.0_dp*cx_i*ux_i(3) - 2.0_dp*zhat) + &
718 qyy_i*(6.0_dp*cy_i*uy_i(3) - 2.0_dp*zhat) + &
719 qzz_i*(6.0_dp*cz_i*uz_i(3) - 2.0_dp*zhat) )
720
721 dudux_i(1) = dudux_i(1) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*xhat)
722 dudux_i(2) = dudux_i(2) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*yhat)
723 dudux_i(3) = dudux_i(3) + pref * sw * ri3 * (qxx_i*6.0_dp*cx_i*zhat)
724
725 duduy_i(1) = duduy_i(1) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*xhat)
726 duduy_i(2) = duduy_i(2) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*yhat)
727 duduy_i(3) = duduy_i(3) + pref * sw * ri3 * (qyy_i*6.0_dp*cy_i*zhat)
728
729 duduz_i(1) = duduz_i(1) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*xhat)
730 duduz_i(2) = duduz_i(2) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*yhat)
731 duduz_i(3) = duduz_i(3) + pref * sw * ri3 * (qzz_i*6.0_dp*cz_i*zhat)
732 endif
733 endif
734
735
736 if (do_pot) then
737 #ifdef IS_MPI
738 pot_row(atom1) = pot_row(atom1) + 0.5d0*epot
739 pot_col(atom2) = pot_col(atom2) + 0.5d0*epot
740 #else
741 pot = pot + epot
742 #endif
743 endif
744
745 #ifdef IS_MPI
746 f_Row(1,atom1) = f_Row(1,atom1) + dudx
747 f_Row(2,atom1) = f_Row(2,atom1) + dudy
748 f_Row(3,atom1) = f_Row(3,atom1) + dudz
749
750 f_Col(1,atom2) = f_Col(1,atom2) - dudx
751 f_Col(2,atom2) = f_Col(2,atom2) - dudy
752 f_Col(3,atom2) = f_Col(3,atom2) - dudz
753
754 if (i_is_Dipole .or. i_is_Quadrupole) then
755 t_Row(1,atom1)=t_Row(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
756 t_Row(2,atom1)=t_Row(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
757 t_Row(3,atom1)=t_Row(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
758 endif
759 if (i_is_Quadrupole) then
760 t_Row(1,atom1)=t_Row(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
761 t_Row(2,atom1)=t_Row(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
762 t_Row(3,atom1)=t_Row(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
763
764 t_Row(1,atom1)=t_Row(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
765 t_Row(2,atom1)=t_Row(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
766 t_Row(3,atom1)=t_Row(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
767 endif
768
769 if (j_is_Dipole .or. j_is_Quadrupole) then
770 t_Col(1,atom2)=t_Col(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
771 t_Col(2,atom2)=t_Col(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
772 t_Col(3,atom2)=t_Col(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
773 endif
774 if (j_is_Quadrupole) then
775 t_Col(1,atom2)=t_Col(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
776 t_Col(2,atom2)=t_Col(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
777 t_Col(3,atom2)=t_Col(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
778
779 t_Col(1,atom2)=t_Col(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
780 t_Col(2,atom2)=t_Col(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
781 t_Col(3,atom2)=t_Col(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
782 endif
783
784 #else
785 f(1,atom1) = f(1,atom1) + dudx
786 f(2,atom1) = f(2,atom1) + dudy
787 f(3,atom1) = f(3,atom1) + dudz
788
789 f(1,atom2) = f(1,atom2) - dudx
790 f(2,atom2) = f(2,atom2) - dudy
791 f(3,atom2) = f(3,atom2) - dudz
792
793 if (i_is_Dipole .or. i_is_Quadrupole) then
794 t(1,atom1)=t(1,atom1) - uz_i(2)*duduz_i(3) + uz_i(3)*duduz_i(2)
795 t(2,atom1)=t(2,atom1) - uz_i(3)*duduz_i(1) + uz_i(1)*duduz_i(3)
796 t(3,atom1)=t(3,atom1) - uz_i(1)*duduz_i(2) + uz_i(2)*duduz_i(1)
797 endif
798 if (i_is_Quadrupole) then
799 t(1,atom1)=t(1,atom1) - ux_i(2)*dudux_i(3) + ux_i(3)*dudux_i(2)
800 t(2,atom1)=t(2,atom1) - ux_i(3)*dudux_i(1) + ux_i(1)*dudux_i(3)
801 t(3,atom1)=t(3,atom1) - ux_i(1)*dudux_i(2) + ux_i(2)*dudux_i(1)
802
803 t(1,atom1)=t(1,atom1) - uy_i(2)*duduy_i(3) + uy_i(3)*duduy_i(2)
804 t(2,atom1)=t(2,atom1) - uy_i(3)*duduy_i(1) + uy_i(1)*duduy_i(3)
805 t(3,atom1)=t(3,atom1) - uy_i(1)*duduy_i(2) + uy_i(2)*duduy_i(1)
806 endif
807
808 if (j_is_Dipole .or. j_is_Quadrupole) then
809 t(1,atom2)=t(1,atom2) - uz_j(2)*duduz_j(3) + uz_j(3)*duduz_j(2)
810 t(2,atom2)=t(2,atom2) - uz_j(3)*duduz_j(1) + uz_j(1)*duduz_j(3)
811 t(3,atom2)=t(3,atom2) - uz_j(1)*duduz_j(2) + uz_j(2)*duduz_j(1)
812 endif
813 if (j_is_Quadrupole) then
814 t(1,atom2)=t(1,atom2) - ux_j(2)*dudux_j(3) + ux_j(3)*dudux_j(2)
815 t(2,atom2)=t(2,atom2) - ux_j(3)*dudux_j(1) + ux_j(1)*dudux_j(3)
816 t(3,atom2)=t(3,atom2) - ux_j(1)*dudux_j(2) + ux_j(2)*dudux_j(1)
817
818 t(1,atom2)=t(1,atom2) - uy_j(2)*duduy_j(3) + uy_j(3)*duduy_j(2)
819 t(2,atom2)=t(2,atom2) - uy_j(3)*duduy_j(1) + uy_j(1)*duduy_j(3)
820 t(3,atom2)=t(3,atom2) - uy_j(1)*duduy_j(2) + uy_j(2)*duduy_j(1)
821 endif
822
823 #endif
824
825 #ifdef IS_MPI
826 id1 = AtomRowToGlobal(atom1)
827 id2 = AtomColToGlobal(atom2)
828 #else
829 id1 = atom1
830 id2 = atom2
831 #endif
832
833 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
834
835 fpair(1) = fpair(1) + dudx
836 fpair(2) = fpair(2) + dudy
837 fpair(3) = fpair(3) + dudz
838
839 endif
840
841 return
842 end subroutine doElectrostaticPair
843
844 !! calculates the switching functions and their derivatives for a given
845 subroutine calc_switch(r, mu, scale, dscale)
846
847 real (kind=dp), intent(in) :: r, mu
848 real (kind=dp), intent(inout) :: scale, dscale
849 real (kind=dp) :: rl, ru, mulow, minRatio, temp, scaleVal
850
851 ! distances must be in angstroms
852 rl = 2.75d0
853 ru = 2.85d0
854 mulow = 3.3856d0 ! 1.84 * 1.84
855 minRatio = mulow / (mu*mu)
856 scaleVal = 1.0d0 - minRatio
857
858 if (r.lt.rl) then
859 scale = minRatio
860 dscale = 0.0d0
861 elseif (r.gt.ru) then
862 scale = 1.0d0
863 dscale = 0.0d0
864 else
865 scale = 1.0d0 - scaleVal*((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2) &
866 / ((ru - rl)**3)
867 dscale = -scaleVal * 6.0d0 * (r-ru)*(r-rl)/((ru - rl)**3)
868 endif
869
870 return
871 end subroutine calc_switch
872
873 subroutine destroyElectrostaticTypes()
874
875 if(allocated(ElectrostaticMap)) deallocate(ElectrostaticMap)
876
877 end subroutine destroyElectrostaticTypes
878
879 end module electrostatic_module