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