ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 1622
Committed: Mon Sep 12 21:49:12 2011 UTC (13 years, 11 months ago) by gezelter
File size: 26934 byte(s)
Log Message:
Adding repulsive power potential

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

Properties

Name Value
svn:keywords Author Id Revision Date