1 |
!! |
2 |
!! Copyright (c) 2007 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. Redistributions of source code must retain the above copyright |
10 |
!! notice, this list of conditions and the following disclaimer. |
11 |
!! |
12 |
!! 2. Redistributions in binary form must reproduce the above copyright |
13 |
!! notice, this list of conditions and the following disclaimer in the |
14 |
!! documentation and/or other materials provided with the |
15 |
!! distribution. |
16 |
!! |
17 |
!! This software is provided "AS IS," without a warranty of any |
18 |
!! kind. All express or implied conditions, representations and |
19 |
!! warranties, including any implied warranty of merchantability, |
20 |
!! fitness for a particular purpose or non-infringement, are hereby |
21 |
!! excluded. The University of Notre Dame and its licensors shall not |
22 |
!! be liable for any damages suffered by licensee as a result of |
23 |
!! using, modifying or distributing the software or its |
24 |
!! derivatives. In no event will the University of Notre Dame or its |
25 |
!! licensors be liable for any lost revenue, profit or data, or for |
26 |
!! direct, indirect, special, consequential, incidental or punitive |
27 |
!! damages, however caused and regardless of the theory of liability, |
28 |
!! arising out of the use of or inability to use software, even if the |
29 |
!! University of Notre Dame has been advised of the possibility of |
30 |
!! such damages. |
31 |
!! |
32 |
!! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
!! research, please cite the appropriate papers when you publish your |
34 |
!! work. Good starting points are: |
35 |
!! |
36 |
!! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
!! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
!! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). |
39 |
!! [4] Vardeman & Gezelter, in progress (2009). |
40 |
!! |
41 |
|
42 |
|
43 |
!! Calculates Metal-Non Metal interactions. |
44 |
!! @author Charles F. Vardeman II |
45 |
!! @version $Id$, $Date$, $Name: not supported by cvs2svn $, $Revision$ |
46 |
|
47 |
|
48 |
module MetalNonMetal |
49 |
use definitions |
50 |
use atype_module |
51 |
use vector_class |
52 |
use simulation |
53 |
use status |
54 |
use fForceOptions |
55 |
use force_globals |
56 |
|
57 |
implicit none |
58 |
PRIVATE |
59 |
#define __FORTRAN90 |
60 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
61 |
#include "UseTheForce/DarkSide/fMnMInteractions.h" |
62 |
|
63 |
logical, save :: useGeometricDistanceMixing = .false. |
64 |
logical, save :: haveInteractionLookup = .false. |
65 |
|
66 |
real(kind=DP), save :: defaultCutoff = 0.0_DP |
67 |
|
68 |
logical, save :: defaultShiftPot = .false. |
69 |
logical, save :: defaultShiftFrc = .false. |
70 |
logical, save :: haveDefaultCutoff = .false. |
71 |
|
72 |
type :: MnMinteraction |
73 |
integer :: metal_atid |
74 |
integer :: nonmetal_atid |
75 |
integer :: interaction_type |
76 |
real(kind=dp) :: R0 |
77 |
real(kind=dp) :: D0 |
78 |
real(kind=dp) :: beta0 |
79 |
real(kind=dp) :: betaH |
80 |
real(kind=dp) :: ca1 |
81 |
real(kind=dp) :: cb1 |
82 |
real(kind=dp) :: sigma |
83 |
real(kind=dp) :: epsilon |
84 |
real(kind=dp) :: rCut = 0.0_dp |
85 |
logical :: rCutWasSet = .false. |
86 |
logical :: shiftedPot = .false. |
87 |
logical :: shiftedFrc = .false. |
88 |
end type MnMinteraction |
89 |
|
90 |
type :: MnMinteractionMap |
91 |
PRIVATE |
92 |
integer :: initialCapacity = 10 |
93 |
integer :: capacityIncrement = 0 |
94 |
integer :: interactionCount = 0 |
95 |
type(MnMinteraction), pointer :: interactions(:) => null() |
96 |
end type MnMinteractionMap |
97 |
|
98 |
type (MnMInteractionMap), pointer :: MnM_Map |
99 |
|
100 |
integer, allocatable, dimension(:,:) :: MnMinteractionLookup |
101 |
|
102 |
public :: setMnMDefaultCutoff |
103 |
public :: addInteraction |
104 |
public :: deleteInteractions |
105 |
public :: MNMtype |
106 |
public :: do_mnm_pair |
107 |
|
108 |
contains |
109 |
|
110 |
|
111 |
subroutine do_mnm_pair(atid1, atid2, D, Rij, R2, Rcut, Sw, vdwMult, & |
112 |
Vpair, Pot, A1, A2, f1, t1, t2) |
113 |
integer, intent(in) :: atid1, atid2 |
114 |
integer :: ljt1, ljt2 |
115 |
real( kind = dp ), intent(in) :: rij, r2, rcut, vdwMult |
116 |
real( kind = dp ) :: pot, sw, vpair |
117 |
real( kind = dp ), intent(inout), dimension(3) :: f1 |
118 |
real (kind=dp), intent(inout), dimension(9) :: A1, A2 |
119 |
real (kind=dp), intent(inout), dimension(3) :: t1, t2 |
120 |
real( kind = dp ), intent(in), dimension(3) :: d |
121 |
|
122 |
integer :: interaction_id |
123 |
integer :: interaction_type |
124 |
|
125 |
if(.not.haveInteractionLookup) then |
126 |
call createInteractionLookup(MnM_MAP) |
127 |
haveInteractionLookup =.true. |
128 |
end if |
129 |
|
130 |
interaction_id = MnMinteractionLookup(atid1, atid2) |
131 |
interaction_type = MnM_Map%interactions(interaction_id)%interaction_type |
132 |
|
133 |
select case (interaction_type) |
134 |
case (MNM_LENNARDJONES) |
135 |
call calc_mnm_lennardjones(D, Rij, R2, Rcut, Sw, & |
136 |
vdwMult, Vpair, Pot, f1, interaction_id) |
137 |
case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE) |
138 |
call calc_mnm_morse(D, Rij, R2, Rcut, Sw, vdwMult, & |
139 |
Vpair, Pot, f1, interaction_id, interaction_type) |
140 |
case(MNM_MAW) |
141 |
call calc_mnm_maw(atid1, atid2, D, Rij, R2, Rcut, Sw, vdwMult, & |
142 |
Vpair, Pot, A1, A2, f1, t1, t2, interaction_id) |
143 |
case default |
144 |
call handleError("MetalNonMetal","Unknown interaction type") |
145 |
end select |
146 |
|
147 |
end subroutine do_mnm_pair |
148 |
|
149 |
subroutine calc_mnm_lennardjones(D, Rij, R2, Rcut, Sw, & |
150 |
vdwMult,Vpair, Pot, f1, interaction_id) |
151 |
|
152 |
real( kind = dp ), intent(in) :: rij, r2, rcut, vdwMult |
153 |
real( kind = dp ) :: pot, sw, vpair |
154 |
real( kind = dp ), intent(inout), dimension(3) :: f1 |
155 |
real( kind = dp ), intent(in), dimension(3) :: d |
156 |
integer, intent(in) :: interaction_id |
157 |
|
158 |
! local Variables |
159 |
real( kind = dp ) :: drdx, drdy, drdz |
160 |
real( kind = dp ) :: fx, fy, fz |
161 |
real( kind = dp ) :: myPot, myPotC, myDeriv, myDerivC, ros, rcos |
162 |
real( kind = dp ) :: pot_temp, dudr |
163 |
real( kind = dp ) :: sigmai |
164 |
real( kind = dp ) :: epsilon |
165 |
logical :: isSoftCore, shiftedPot, shiftedFrc |
166 |
integer :: id1, id2, localError |
167 |
|
168 |
sigmai = 1.0_dp / MnM_Map%interactions(interaction_id)%sigma |
169 |
epsilon = MnM_Map%interactions(interaction_id)%epsilon |
170 |
shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot |
171 |
shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc |
172 |
|
173 |
ros = rij * sigmai |
174 |
|
175 |
call getLJfunc(ros, myPot, myDeriv) |
176 |
|
177 |
if (shiftedPot) then |
178 |
rcos = rcut * sigmai |
179 |
call getLJfunc(rcos, myPotC, myDerivC) |
180 |
myDerivC = 0.0_dp |
181 |
elseif (shiftedFrc) then |
182 |
rcos = rcut * sigmai |
183 |
call getLJfunc(rcos, myPotC, myDerivC) |
184 |
myPotC = myPotC + myDerivC * (rij - rcut) * sigmai |
185 |
else |
186 |
myPotC = 0.0_dp |
187 |
myDerivC = 0.0_dp |
188 |
endif |
189 |
|
190 |
pot_temp = vdwMult * epsilon * (myPot - myPotC) |
191 |
vpair = vpair + pot_temp |
192 |
dudr = sw * vdwMult * epsilon * (myDeriv - myDerivC) * sigmai |
193 |
|
194 |
drdx = d(1) / rij |
195 |
drdy = d(2) / rij |
196 |
drdz = d(3) / rij |
197 |
|
198 |
fx = dudr * drdx |
199 |
fy = dudr * drdy |
200 |
fz = dudr * drdz |
201 |
|
202 |
pot = pot + sw*pot_temp |
203 |
f1(1) = f1(1) + fx |
204 |
f1(2) = f1(2) + fy |
205 |
f1(3) = f1(3) + fz |
206 |
|
207 |
return |
208 |
end subroutine calc_mnm_lennardjones |
209 |
|
210 |
subroutine calc_mnm_morse(D, Rij, R2, Rcut, Sw, vdwMult, & |
211 |
Vpair, Pot, f1, interaction_id, interaction_type) |
212 |
real( kind = dp ), intent(in) :: rij, r2, rcut, vdwMult |
213 |
real( kind = dp ) :: pot, sw, vpair |
214 |
real( kind = dp ), intent(inout), dimension(3) :: f1 |
215 |
real( kind = dp ), intent(in), dimension(3) :: d |
216 |
integer, intent(in) :: interaction_id, interaction_type |
217 |
logical :: shiftedPot, shiftedFrc |
218 |
|
219 |
! Local Variables |
220 |
real(kind=dp) :: Beta0 |
221 |
real(kind=dp) :: R0 |
222 |
real(kind=dp) :: D0 |
223 |
real(kind=dp) :: expt |
224 |
real(kind=dp) :: expt2 |
225 |
real(kind=dp) :: expfnc |
226 |
real(kind=dp) :: expfnc2 |
227 |
real(kind=dp) :: D_expt |
228 |
real(kind=dp) :: D_expt2 |
229 |
real(kind=dp) :: rcos |
230 |
real(kind=dp) :: exptC |
231 |
real(kind=dp) :: expt2C |
232 |
real(kind=dp) :: expfncC |
233 |
real(kind=dp) :: expfnc2C |
234 |
real(kind=dp) :: D_expt2C |
235 |
real(kind=dp) :: D_exptC |
236 |
|
237 |
real(kind=dp) :: myPot |
238 |
real(kind=dp) :: myPotC |
239 |
real(kind=dp) :: myDeriv |
240 |
real(kind=dp) :: myDerivC |
241 |
real(kind=dp) :: pot_temp |
242 |
real(kind=dp) :: fx,fy,fz |
243 |
real(kind=dp) :: drdx,drdy,drdz |
244 |
real(kind=dp) :: dudr |
245 |
integer :: id1,id2 |
246 |
|
247 |
|
248 |
D0 = MnM_Map%interactions(interaction_id)%D0 |
249 |
R0 = MnM_Map%interactions(interaction_id)%r0 |
250 |
Beta0 = MnM_Map%interactions(interaction_id)%Beta0 |
251 |
shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot |
252 |
shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc |
253 |
|
254 |
! V(r) = D_e exp(-a(r-re)(exp(-a(r-re))-2) |
255 |
|
256 |
expt = -beta0*(rij - R0) |
257 |
expfnc = exp(expt) |
258 |
expfnc2 = expfnc*expfnc |
259 |
|
260 |
if (shiftedPot .or. shiftedFrc) then |
261 |
exptC = -beta0*(rcut - R0) |
262 |
expfncC = exp(exptC) |
263 |
expfnc2C = expfncC*expfncC |
264 |
endif |
265 |
|
266 |
select case (interaction_type) |
267 |
case (MNM_SHIFTEDMORSE) |
268 |
|
269 |
myPot = D0 * (expfnc2 - 2.0_dp * expfnc) |
270 |
myDeriv = 2.0*D0*beta0*(expfnc - expfnc2) |
271 |
|
272 |
if (shiftedPot) then |
273 |
myPotC = D0 * (expfnc2C - 2.0_dp * expfncC) |
274 |
myDerivC = 0.0_dp |
275 |
elseif (shiftedFrc) then |
276 |
myPotC = D0 * (expfnc2C - 2.0_dp * expfncC) |
277 |
myDerivC = 2.0*D0*beta0*(expfnc2C - expfnc2C) |
278 |
myPotC = myPotC + myDerivC * (rij - rcut) |
279 |
else |
280 |
myPotC = 0.0_dp |
281 |
myDerivC = 0.0_dp |
282 |
endif |
283 |
|
284 |
case (MNM_REPULSIVEMORSE) |
285 |
|
286 |
myPot = D0 * expfnc2 |
287 |
myDeriv = -2.0_dp * D0 * beta0 * expfnc2 |
288 |
|
289 |
if (shiftedPot) then |
290 |
myPotC = D0 * expfnc2C |
291 |
myDerivC = 0.0_dp |
292 |
elseif (shiftedFrc) then |
293 |
myPotC = D0 * expfnc2C |
294 |
myDerivC = -2.0_dp * D0* beta0 * expfnc2C |
295 |
myPotC = myPotC + myDerivC * (rij - rcut) |
296 |
else |
297 |
myPotC = 0.0_dp |
298 |
myDerivC = 0.0_dp |
299 |
endif |
300 |
end select |
301 |
|
302 |
pot_temp = vdwMult * (myPot - myPotC) |
303 |
vpair = vpair + pot_temp |
304 |
dudr = sw * vdwMult * (myDeriv - myDerivC) |
305 |
|
306 |
drdx = d(1) / rij |
307 |
drdy = d(2) / rij |
308 |
drdz = d(3) / rij |
309 |
|
310 |
fx = dudr * drdx |
311 |
fy = dudr * drdy |
312 |
fz = dudr * drdz |
313 |
|
314 |
pot = pot + sw*pot_temp |
315 |
|
316 |
f1(1) = f1(1) + fx |
317 |
f1(2) = f1(2) + fy |
318 |
f1(3) = f1(3) + fz |
319 |
|
320 |
return |
321 |
end subroutine calc_mnm_morse |
322 |
|
323 |
subroutine calc_mnm_maw(atid1, atid2, D, Rij, R2, Rcut, Sw, vdwMult, & |
324 |
Vpair, Pot, A1, A2, f1, t1, t2, interaction_id) |
325 |
real( kind = dp ), intent(in) :: rij, r2, rcut, vdwMult |
326 |
real( kind = dp ) :: pot, sw, vpair |
327 |
real( kind = dp ), intent(inout), dimension(3) :: f1 |
328 |
real (kind=dp),intent(in), dimension(9) :: A1, A2 |
329 |
real (kind=dp),intent(inout), dimension(3) :: t1, t2 |
330 |
|
331 |
real( kind = dp ), intent(in), dimension(3) :: d |
332 |
|
333 |
integer, intent(in) :: interaction_id |
334 |
|
335 |
real(kind = dp) :: D0, R0, beta0, ca1, cb1, pot_temp |
336 |
real(kind = dp) :: expt0, expfnc0, expfnc02 |
337 |
real(kind = dp) :: exptH, expfncH, expfncH2 |
338 |
real(kind = dp) :: x, y, z, x2, y2, z2, r3, r4 |
339 |
real(kind = dp) :: drdx, drdy, drdz |
340 |
real(kind = dp) :: dvdx, dvdy, dvdz |
341 |
real(kind = dp) :: Vang, dVangdx, dVangdy, dVangdz |
342 |
real(kind = dp) :: dVangdux, dVangduy, dVangduz |
343 |
real(kind = dp) :: dVmorsedr |
344 |
real(kind = dp) :: Vmorse, dVmorsedx, dVmorsedy, dVmorsedz |
345 |
real(kind = dp) :: ta1, b1, s |
346 |
real(kind = dp) :: da1dx, da1dy, da1dz, da1dux, da1duy, da1duz |
347 |
real(kind = dp) :: db1dx, db1dy, db1dz, db1dux, db1duy, db1duz |
348 |
real(kind = dp) :: fx, fy, fz, tx, ty, tz, fxl, fyl, fzl |
349 |
! real(kind = dp), parameter :: st = sqrt(3.0_dp) |
350 |
real(kind = dp), parameter :: st = 1.732050807568877 |
351 |
integer :: atid1, atid2, id1, id2 |
352 |
logical :: shiftedPot, shiftedFrc |
353 |
|
354 |
if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then |
355 |
! rotate the inter-particle separation into the two different |
356 |
! body-fixed coordinate systems: |
357 |
|
358 |
x = A1(1)*d(1) + A1(2)*d(2) + A1(3)*d(3) |
359 |
y = A1(4)*d(1) + A1(5)*d(2) + A1(6)*d(3) |
360 |
z = A1(7)*d(1) + A1(8)*d(2) + A1(9)*d(3) |
361 |
else |
362 |
! negative sign because this is the vector from j to i: |
363 |
|
364 |
x = -(A2(1)*d(1) + A2(2)*d(2) + A2(3)*d(3)) |
365 |
y = -(A2(4)*d(1) + A2(5)*d(2) + A2(6)*d(3)) |
366 |
z = -(A2(7)*d(1) + A2(8)*d(2) + A2(9)*d(3)) |
367 |
endif |
368 |
|
369 |
D0 = MnM_Map%interactions(interaction_id)%D0 |
370 |
R0 = MnM_Map%interactions(interaction_id)%r0 |
371 |
beta0 = MnM_Map%interactions(interaction_id)%beta0 |
372 |
ca1 = MnM_Map%interactions(interaction_id)%ca1 |
373 |
cb1 = MnM_Map%interactions(interaction_id)%cb1 |
374 |
|
375 |
shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot |
376 |
shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc |
377 |
|
378 |
expt0 = -beta0*(rij - R0) |
379 |
expfnc0 = exp(expt0) |
380 |
expfnc02 = expfnc0*expfnc0 |
381 |
|
382 |
!!$ if (shiftedPot .or. shiftedFrc) then |
383 |
!!$ exptC0 = -beta0*(rcut - R0) |
384 |
!!$ expfncC0 = exp(exptC0) |
385 |
!!$ expfncC02 = expfncC0*expfncC0 |
386 |
!!$ exptCH = -betaH*(rcut - R0) |
387 |
!!$ expfncCH = exp(exptCH) |
388 |
!!$ expfncCH2 = expfncCH*expfncCH |
389 |
!!$ endif |
390 |
|
391 |
drdx = x / rij |
392 |
drdy = y / rij |
393 |
drdz = z / rij |
394 |
|
395 |
x2 = x*x |
396 |
y2 = y*y |
397 |
z2 = z*z |
398 |
r3 = r2*rij |
399 |
r4 = r2*r2 |
400 |
|
401 |
Vmorse = D0 * (expfnc02 - 2.0_dp * expfnc0) |
402 |
|
403 |
! angular modulation of morse part of potential to approximate |
404 |
! the squares of the two HOMO lone pair orbitals in water: |
405 |
!********************** old form************************* |
406 |
! s = 1 / (4 pi) |
407 |
! ta1 = (s - pz)^2 = (1 - sqrt(3)*cos(theta))^2 / (4 pi) |
408 |
! b1 = px^2 = 3 * (sin(theta)*cos(phi))^2 / (4 pi) |
409 |
!********************** old form************************* |
410 |
! we'll leave out the 4 pi for now |
411 |
|
412 |
! new functional form just using the p orbitals. |
413 |
! Vmorse(r)*[a*p_x + b p_z + (1-a-b)] |
414 |
! which is |
415 |
! Vmorse(r)*[a sin^2(theta) cos^2(phi) + b cos(theta) + (1-a-b)] |
416 |
! Vmorse(r)*[a*x2/r2 + b*z/r + (1-a-b)] |
417 |
|
418 |
|
419 |
|
420 |
s = 1.0_dp |
421 |
! ta1 = (1.0_dp - st * z / rij )**2 |
422 |
! b1 = 3.0_dp * x2 / r2 |
423 |
|
424 |
! Vang = s + ca1 * ta1 + cb1 * b1 |
425 |
|
426 |
Vang = ca1 * x2/r2 + cb1 * z/rij + (0.8_dp-ca1/3.0_dp) |
427 |
|
428 |
pot_temp = vdwMult * Vmorse*Vang |
429 |
|
430 |
vpair = vpair + pot_temp |
431 |
pot = pot + pot_temp*sw |
432 |
|
433 |
dVmorsedr = 2.0_dp*D0*beta0*(expfnc0 - expfnc02) |
434 |
|
435 |
dVmorsedx = dVmorsedr * drdx |
436 |
dVmorsedy = dVmorsedr * drdy |
437 |
dVmorsedz = dVmorsedr * drdz |
438 |
|
439 |
!da1dx = 2.0_dp * st * x * z / r3 - 6.0_dp * x * z2 / r4 |
440 |
!da1dy = 2.0_dp * st * y * z / r3 - 6.0_dp * y * z2 / r4 |
441 |
!da1dz = 2.0_dp * st * (x2 + y2) * (st * z - rij ) / r4 |
442 |
|
443 |
!db1dx = 6.0_dp * x * (1.0_dp - x2 / r2) / r2 |
444 |
!db1dy = -6.0_dp * x2 * y / r4 |
445 |
!db1dz = -6.0_dp * x2 * z / r4 |
446 |
|
447 |
!dVangdx = ca1 * da1dx + cb1 * db1dx |
448 |
!dVangdy = ca1 * da1dy + cb1 * db1dy |
449 |
!dVangdz = ca1 * da1dz + cb1 * db1dz |
450 |
|
451 |
dVangdx = -2.0*ca1*x2*x/r4 + 2.0*ca1*x/r2 - cb1*x*z/r3 |
452 |
dVangdy = -2.0*ca1*x2*y/r4 - cb1*z*y/r3 |
453 |
dVangdz = -2.0*ca1*x2*z/r4 + cb1/rij - cb1*z2 /r3 |
454 |
|
455 |
! chain rule to put these back on x, y, z |
456 |
dvdx = Vang * dVmorsedx + Vmorse * dVangdx |
457 |
dvdy = Vang * dVmorsedy + Vmorse * dVangdy |
458 |
dvdz = Vang * dVmorsedz + Vmorse * dVangdz |
459 |
|
460 |
! Torques for Vang using method of Price: |
461 |
! S. L. Price, A. J. Stone, and M. Alderton, Mol. Phys. 52, 987 (1984). |
462 |
|
463 |
!da1dux = 6.0_dp * y * z / r2 - 2.0_dp * st * y / rij |
464 |
!da1duy = -6.0_dp * x * z / r2 + 2.0_dp * st * x / rij |
465 |
!da1duz = 0.0_dp |
466 |
|
467 |
!db1dux = 0.0_dp |
468 |
!db1duy = 6.0_dp * x * z / r2 |
469 |
!db1duz = -6.0_dp * y * x / r2 |
470 |
|
471 |
!dVangdux = ca1 * da1dux + cb1 * db1dux |
472 |
!dVangduy = ca1 * da1duy + cb1 * db1duy |
473 |
!dVangduz = ca1 * da1duz + cb1 * db1duz |
474 |
|
475 |
dVangdux = cb1*y/rij |
476 |
dVangduy = 2.0*ca1*x*z/r2 - cb1*x/rij |
477 |
dVangduz = -2.0*ca1*y*x/r2 |
478 |
|
479 |
! do the torques first since they are easy: |
480 |
! remember that these are still in the body fixed axes |
481 |
|
482 |
tx = vdwMult * Vmorse * dVangdux * sw |
483 |
ty = vdwMult * Vmorse * dVangduy * sw |
484 |
tz = vdwMult * Vmorse * dVangduz * sw |
485 |
|
486 |
! go back to lab frame using transpose of rotation matrix: |
487 |
|
488 |
if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then |
489 |
t1(1) = t1(1) + a1(1)*tx + a1(4)*ty + a1(7)*tz |
490 |
t1(2) = t1(2) + a1(2)*tx + a1(5)*ty + a1(8)*tz |
491 |
t1(3) = t1(3) + a1(3)*tx + a1(6)*ty + a1(9)*tz |
492 |
else |
493 |
t2(1) = t2(1) + a2(1)*tx + a2(4)*ty + a2(7)*tz |
494 |
t2(2) = t2(2) + a2(2)*tx + a2(5)*ty + a2(8)*tz |
495 |
t2(3) = t2(3) + a2(3)*tx + a2(6)*ty + a2(9)*tz |
496 |
endif |
497 |
|
498 |
! Now, on to the forces (still in body frame of water) |
499 |
|
500 |
fx = vdwMult * dvdx * sw |
501 |
fy = vdwMult * dvdy * sw |
502 |
fz = vdwMult * dvdz * sw |
503 |
|
504 |
! rotate the terms back into the lab frame: |
505 |
|
506 |
if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then |
507 |
fxl = a1(1)*fx + a1(4)*fy + a1(7)*fz |
508 |
fyl = a1(2)*fx + a1(5)*fy + a1(8)*fz |
509 |
fzl = a1(3)*fx + a1(6)*fy + a1(9)*fz |
510 |
else |
511 |
! negative sign because this is the vector from j to i: |
512 |
fxl = -(a2(1)*fx + a2(4)*fy + a2(7)*fz) |
513 |
fyl = -(a2(2)*fx + a2(5)*fy + a2(8)*fz) |
514 |
fzl = -(a2(3)*fx + a2(6)*fy + a2(9)*fz) |
515 |
endif |
516 |
|
517 |
f1(1) = f1(1) + fxl |
518 |
f1(2) = f1(2) + fyl |
519 |
f1(3) = f1(3) + fzl |
520 |
|
521 |
return |
522 |
end subroutine calc_mnm_maw |
523 |
|
524 |
|
525 |
subroutine setMnMDefaultCutoff(thisRcut, shiftedPot, shiftedFrc) |
526 |
real(kind=dp), intent(in) :: thisRcut |
527 |
logical, intent(in) :: shiftedPot |
528 |
logical, intent(in) :: shiftedFrc |
529 |
integer i, nInteractions |
530 |
defaultCutoff = thisRcut |
531 |
defaultShiftPot = shiftedPot |
532 |
defaultShiftFrc = shiftedFrc |
533 |
|
534 |
if (associated(MnM_Map)) then |
535 |
if(MnM_Map%interactionCount /= 0) then |
536 |
nInteractions = MnM_Map%interactionCount |
537 |
|
538 |
do i = 1, nInteractions |
539 |
MnM_Map%interactions(i)%shiftedPot = shiftedPot |
540 |
MnM_Map%interactions(i)%shiftedFrc = shiftedFrc |
541 |
MnM_Map%interactions(i)%rCut = thisRcut |
542 |
MnM_Map%interactions(i)%rCutWasSet = .true. |
543 |
enddo |
544 |
end if |
545 |
end if |
546 |
|
547 |
end subroutine setMnMDefaultCutoff |
548 |
|
549 |
subroutine copyAllData(v1, v2) |
550 |
type(MnMinteractionMap), pointer :: v1 |
551 |
type(MnMinteractionMap), pointer :: v2 |
552 |
integer :: i, j |
553 |
|
554 |
do i = 1, v1%interactionCount |
555 |
v2%interactions(i) = v1%interactions(i) |
556 |
enddo |
557 |
|
558 |
v2%interactionCount = v1%interactionCount |
559 |
return |
560 |
end subroutine copyAllData |
561 |
|
562 |
subroutine addInteraction(myInteraction) |
563 |
type(MNMtype) :: myInteraction |
564 |
type(MnMinteraction) :: nt |
565 |
integer :: id |
566 |
|
567 |
nt%interaction_type = myInteraction%MNMInteractionType |
568 |
nt%metal_atid = & |
569 |
getFirstMatchingElement(atypes, "c_ident", myInteraction%metal_atid) |
570 |
nt%nonmetal_atid = & |
571 |
getFirstMatchingElement(atypes, "c_ident", myInteraction%nonmetal_atid) |
572 |
|
573 |
select case (nt%interaction_type) |
574 |
case (MNM_LENNARDJONES) |
575 |
nt%sigma = myInteraction%sigma |
576 |
nt%epsilon = myInteraction%epsilon |
577 |
case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE) |
578 |
nt%R0 = myInteraction%R0 |
579 |
nt%D0 = myInteraction%D0 |
580 |
nt%beta0 = myInteraction%beta0 |
581 |
case(MNM_MAW) |
582 |
nt%R0 = myInteraction%R0 |
583 |
nt%D0 = myInteraction%D0 |
584 |
nt%beta0 = myInteraction%beta0 |
585 |
nt%ca1 = myInteraction%ca1 |
586 |
nt%cb1 = myInteraction%cb1 |
587 |
case default |
588 |
call handleError("MNM", "Unknown Interaction type") |
589 |
end select |
590 |
|
591 |
if (.not. associated(MnM_Map)) then |
592 |
call ensureCapacityHelper(MnM_Map, 1) |
593 |
else |
594 |
call ensureCapacityHelper(MnM_Map, MnM_Map%interactionCount + 1) |
595 |
end if |
596 |
|
597 |
MnM_Map%interactionCount = MnM_Map%interactionCount + 1 |
598 |
id = MnM_Map%interactionCount |
599 |
MnM_Map%interactions(id) = nt |
600 |
end subroutine addInteraction |
601 |
|
602 |
subroutine ensureCapacityHelper(this, minCapacity) |
603 |
type(MnMinteractionMap), pointer :: this, that |
604 |
integer, intent(in) :: minCapacity |
605 |
integer :: oldCapacity |
606 |
integer :: newCapacity |
607 |
logical :: resizeFlag |
608 |
|
609 |
resizeFlag = .false. |
610 |
|
611 |
! first time: allocate a new vector with default size |
612 |
|
613 |
if (.not. associated(this)) then |
614 |
this => MnMinitialize(minCapacity, 0) |
615 |
endif |
616 |
|
617 |
oldCapacity = size(this%interactions) |
618 |
|
619 |
if (minCapacity > oldCapacity) then |
620 |
if (this%capacityIncrement .gt. 0) then |
621 |
newCapacity = oldCapacity + this%capacityIncrement |
622 |
else |
623 |
newCapacity = oldCapacity * 2 |
624 |
endif |
625 |
if (newCapacity .lt. minCapacity) then |
626 |
newCapacity = minCapacity |
627 |
endif |
628 |
resizeFlag = .true. |
629 |
else |
630 |
newCapacity = oldCapacity |
631 |
endif |
632 |
|
633 |
if (resizeFlag) then |
634 |
that => MnMinitialize(newCapacity, this%capacityIncrement) |
635 |
call copyAllData(this, that) |
636 |
this => MnMdestroy(this) |
637 |
this => that |
638 |
endif |
639 |
end subroutine ensureCapacityHelper |
640 |
|
641 |
function MnMinitialize(cap, capinc) result(this) |
642 |
integer, intent(in) :: cap, capinc |
643 |
integer :: error |
644 |
type(MnMinteractionMap), pointer :: this |
645 |
|
646 |
nullify(this) |
647 |
|
648 |
if (cap < 0) then |
649 |
write(*,*) 'Bogus Capacity:', cap |
650 |
return |
651 |
endif |
652 |
allocate(this,stat=error) |
653 |
if ( error /= 0 ) then |
654 |
write(*,*) 'Could not allocate MnMinteractionMap!' |
655 |
return |
656 |
end if |
657 |
|
658 |
this%initialCapacity = cap |
659 |
this%capacityIncrement = capinc |
660 |
|
661 |
allocate(this%interactions(this%initialCapacity), stat=error) |
662 |
if(error /= 0) write(*,*) 'Could not allocate MnMinteraction!' |
663 |
|
664 |
end function MnMinitialize |
665 |
|
666 |
subroutine createInteractionLookup(this) |
667 |
type (MnMInteractionMap),pointer :: this |
668 |
integer :: biggestAtid, i, metal_atid, nonmetal_atid, error |
669 |
|
670 |
biggestAtid=-1 |
671 |
do i = 1, this%interactionCount |
672 |
metal_atid = this%interactions(i)%metal_atid |
673 |
nonmetal_atid = this%interactions(i)%nonmetal_atid |
674 |
|
675 |
if (metal_atid .gt. biggestAtid) biggestAtid = metal_atid |
676 |
if (nonmetal_atid .gt. biggestAtid) biggestAtid = nonmetal_atid |
677 |
enddo |
678 |
|
679 |
allocate(MnMinteractionLookup(biggestAtid,biggestAtid), stat=error) |
680 |
if (error /= 0) write(*,*) 'Could not allocate MnMinteractionLookup' |
681 |
|
682 |
do i = 1, this%interactionCount |
683 |
metal_atid = this%interactions(i)%metal_atid |
684 |
nonmetal_atid = this%interactions(i)%nonmetal_atid |
685 |
|
686 |
MnMinteractionLookup(metal_atid, nonmetal_atid) = i |
687 |
MnMinteractionLookup(nonmetal_atid, metal_atid) = i |
688 |
enddo |
689 |
end subroutine createInteractionLookup |
690 |
|
691 |
function MnMdestroy(this) result(null_this) |
692 |
logical :: done |
693 |
type(MnMinteractionMap), pointer :: this |
694 |
type(MnMinteractionMap), pointer :: null_this |
695 |
|
696 |
if (.not. associated(this)) then |
697 |
null_this => null() |
698 |
return |
699 |
end if |
700 |
|
701 |
!! Walk down the list and deallocate each of the map's components |
702 |
if(associated(this%interactions)) then |
703 |
deallocate(this%interactions) |
704 |
this%interactions=>null() |
705 |
endif |
706 |
deallocate(this) |
707 |
this => null() |
708 |
null_this => null() |
709 |
end function MnMdestroy |
710 |
|
711 |
subroutine deleteInteractions() |
712 |
MnM_Map => MnMdestroy(MnM_Map) |
713 |
return |
714 |
end subroutine deleteInteractions |
715 |
|
716 |
subroutine getLJfunc(r, myPot, myDeriv) |
717 |
|
718 |
real(kind=dp), intent(in) :: r |
719 |
real(kind=dp), intent(inout) :: myPot, myDeriv |
720 |
real(kind=dp) :: ri, ri2, ri6, ri7, ri12, ri13 |
721 |
real(kind=dp) :: a, b, c, d, dx |
722 |
integer :: j |
723 |
|
724 |
ri = 1.0_DP / r |
725 |
ri2 = ri*ri |
726 |
ri6 = ri2*ri2*ri2 |
727 |
ri7 = ri6*ri |
728 |
ri12 = ri6*ri6 |
729 |
ri13 = ri12*ri |
730 |
|
731 |
myPot = 4.0_DP * (ri12 - ri6) |
732 |
myDeriv = 24.0_DP * (ri7 - 2.0_DP * ri13) |
733 |
|
734 |
return |
735 |
end subroutine getLJfunc |
736 |
|
737 |
subroutine getSoftFunc(r, myPot, myDeriv) |
738 |
|
739 |
real(kind=dp), intent(in) :: r |
740 |
real(kind=dp), intent(inout) :: myPot, myDeriv |
741 |
real(kind=dp) :: ri, ri2, ri6, ri7 |
742 |
real(kind=dp) :: a, b, c, d, dx |
743 |
integer :: j |
744 |
|
745 |
ri = 1.0_DP / r |
746 |
ri2 = ri*ri |
747 |
ri6 = ri2*ri2*ri2 |
748 |
ri7 = ri6*ri |
749 |
myPot = 4.0_DP * (ri6) |
750 |
myDeriv = - 24.0_DP * ri7 |
751 |
|
752 |
return |
753 |
end subroutine getSoftFunc |
754 |
end module MetalNonMetal |