ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 1465
Committed: Fri Jul 9 23:08:25 2010 UTC (15 years ago) by chuckv
File size: 24230 byte(s)
Log Message:
Creating busticated version of OpenMD

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date