ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 1245
Committed: Tue May 27 16:39:06 2008 UTC (17 years, 1 month ago) by chuckv
Original Path: trunk/src/UseTheForce/DarkSide/MetalNonMetal.F90
File size: 30069 byte(s)
Log Message:
Checking in changes for Hefland moment calculations

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     !! 1. Acknowledgement of the program authors must be made in any
10     !! publication of scientific results based in part on use of the
11     !! program. An acceptable form of acknowledgement is citation of
12     !! the article in which the program was described (Matthew
13     !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14     !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     !! Parallel Simulation Engine for Molecular Dynamics,"
16     !! J. Comput. Chem. 26, pp. 252-271 (2005))
17     !!
18     !! 2. Redistributions of source code must retain the above copyright
19     !! notice, this list of conditions and the following disclaimer.
20     !!
21     !! 3. Redistributions in binary form must reproduce the above copyright
22     !! notice, this list of conditions and the following disclaimer in the
23     !! documentation and/or other materials provided with the
24     !! distribution.
25     !!
26     !! This software is provided "AS IS," without a warranty of any
27     !! kind. All express or implied conditions, representations and
28     !! warranties, including any implied warranty of merchantability,
29     !! fitness for a particular purpose or non-infringement, are hereby
30     !! excluded. The University of Notre Dame and its licensors shall not
31     !! be liable for any damages suffered by licensee as a result of
32     !! using, modifying or distributing the software or its
33     !! derivatives. In no event will the University of Notre Dame or its
34     !! licensors be liable for any lost revenue, profit or data, or for
35     !! direct, indirect, special, consequential, incidental or punitive
36     !! damages, however caused and regardless of the theory of liability,
37     !! arising out of the use of or inability to use software, even if the
38     !! University of Notre Dame has been advised of the possibility of
39     !! such damages.
40     !!
41    
42    
43     !! Calculates Metal-Non Metal interactions.
44     !! @author Charles F. Vardeman II
45 chuckv 1245 !! @version $Id: MetalNonMetal.F90,v 1.13 2008-05-27 16:39:05 chuckv Exp $, $Date: 2008-05-27 16:39:05 $, $Name: not supported by cvs2svn $, $Revision: 1.13 $
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     #ifdef IS_MPI
56     use mpiSimulation
57     #endif
58     use force_globals
59    
60     implicit none
61     PRIVATE
62     #define __FORTRAN90
63     #include "UseTheForce/DarkSide/fInteractionMap.h"
64     #include "UseTheForce/DarkSide/fMnMInteractions.h"
65    
66     logical, save :: useGeometricDistanceMixing = .false.
67     logical, save :: haveInteractionLookup = .false.
68    
69     real(kind=DP), save :: defaultCutoff = 0.0_DP
70 chuckv 1228
71 chuckv 1168 logical, save :: defaultShiftPot = .false.
72     logical, save :: defaultShiftFrc = .false.
73     logical, save :: haveDefaultCutoff = .false.
74    
75     type :: MnMinteraction
76     integer :: metal_atid
77     integer :: nonmetal_atid
78     integer :: interaction_type
79     real(kind=dp) :: R0
80     real(kind=dp) :: D0
81     real(kind=dp) :: beta0
82     real(kind=dp) :: betaH
83 gezelter 1238 real(kind=dp) :: ca1
84     real(kind=dp) :: cb1
85 chuckv 1168 real(kind=dp) :: sigma
86     real(kind=dp) :: epsilon
87     real(kind=dp) :: rCut = 0.0_dp
88     logical :: rCutWasSet = .false.
89     logical :: shiftedPot = .false.
90     logical :: shiftedFrc = .false.
91     end type MnMinteraction
92    
93     type :: MnMinteractionMap
94     PRIVATE
95     integer :: initialCapacity = 10
96     integer :: capacityIncrement = 0
97     integer :: interactionCount = 0
98     type(MnMinteraction), pointer :: interactions(:) => null()
99     end type MnMinteractionMap
100    
101 xsun 1178 type (MnMInteractionMap), pointer :: MnM_Map
102 chuckv 1168
103     integer, allocatable, dimension(:,:) :: MnMinteractionLookup
104    
105     public :: setMnMDefaultCutoff
106     public :: addInteraction
107     public :: deleteInteractions
108     public :: MNMtype
109 chuckv 1173 public :: do_mnm_pair
110 chuckv 1168
111     contains
112    
113    
114 chuckv 1173 subroutine do_mnm_pair(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
115     Pot, A, F,t, Do_pot)
116 gezelter 1174 integer, intent(in) :: atom1, atom2
117 chuckv 1173 integer :: atid1, atid2, ljt1, ljt2
118 gezelter 1174 real( kind = dp ), intent(in) :: rij, r2, rcut
119     real( kind = dp ) :: pot, sw, vpair
120 chuckv 1173 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
121 gezelter 1174 real (kind=dp), intent(in), dimension(9,nLocal) :: A
122     real (kind=dp), intent(inout), dimension(3,nLocal) :: t
123     real( kind = dp ), intent(in), dimension(3) :: d
124 chuckv 1173 real( kind = dp ), intent(inout), dimension(3) :: fpair
125 gezelter 1174 logical, intent(in) :: do_pot
126 chuckv 1168
127 chuckv 1173 integer :: interaction_id
128     integer :: interaction_type
129    
130     #ifdef IS_MPI
131     atid1 = atid_Row(atom1)
132     atid2 = atid_Col(atom2)
133     #else
134     atid1 = atid(atom1)
135     atid2 = atid(atom2)
136     #endif
137    
138 chuckv 1175 if(.not.haveInteractionLookup) then
139     call createInteractionLookup(MnM_MAP)
140     haveInteractionLookup =.true.
141     end if
142    
143 chuckv 1173 interaction_id = MnMinteractionLookup(atid1, atid2)
144     interaction_type = MnM_Map%interactions(interaction_id)%interaction_type
145 chuckv 1229
146 chuckv 1173
147 chuckv 1229
148 chuckv 1173 select case (interaction_type)
149     case (MNM_LENNARDJONES)
150 gezelter 1174 call calc_mnm_lennardjones(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, &
151     Fpair, Pot, F, Do_pot, interaction_id)
152 chuckv 1173 case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE)
153     call calc_mnm_morse(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
154 gezelter 1174 Pot, F, Do_pot, interaction_id, interaction_type)
155 chuckv 1173 case(MNM_MAW)
156     call calc_mnm_maw(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
157 gezelter 1174 Pot,A, F,t, Do_pot, interaction_id)
158 chuckv 1229 case default
159     call handleError("MetalNonMetal","Unknown interaction type")
160 chuckv 1173 end select
161    
162     end subroutine do_mnm_pair
163    
164 gezelter 1174 subroutine calc_mnm_lennardjones(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, &
165     Fpair, Pot, F, Do_pot, interaction_id)
166 chuckv 1173
167 gezelter 1174 integer, intent(in) :: atom1, atom2
168     real( kind = dp ), intent(in) :: rij, r2, rcut
169     real( kind = dp ) :: pot, sw, vpair
170 chuckv 1173 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
171 gezelter 1174 real( kind = dp ), intent(in), dimension(3) :: d
172 chuckv 1173 real( kind = dp ), intent(inout), dimension(3) :: fpair
173 gezelter 1174 logical, intent(in) :: do_pot
174 chuckv 1173 integer, intent(in) :: interaction_id
175    
176     ! local Variables
177     real( kind = dp ) :: drdx, drdy, drdz
178     real( kind = dp ) :: fx, fy, fz
179     real( kind = dp ) :: myPot, myPotC, myDeriv, myDerivC, ros, rcos
180     real( kind = dp ) :: pot_temp, dudr
181     real( kind = dp ) :: sigmai
182     real( kind = dp ) :: epsilon
183     logical :: isSoftCore, shiftedPot, shiftedFrc
184     integer :: id1, id2, localError
185    
186 gezelter 1174 sigmai = 1.0_dp / MnM_Map%interactions(interaction_id)%sigma
187 chuckv 1173 epsilon = MnM_Map%interactions(interaction_id)%epsilon
188     shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot
189     shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc
190    
191     ros = rij * sigmai
192    
193     call getLJfunc(ros, myPot, myDeriv)
194    
195     if (shiftedPot) then
196     rcos = rcut * sigmai
197     call getLJfunc(rcos, myPotC, myDerivC)
198     myDerivC = 0.0_dp
199     elseif (shiftedFrc) then
200     rcos = rcut * sigmai
201     call getLJfunc(rcos, myPotC, myDerivC)
202     myPotC = myPotC + myDerivC * (rij - rcut) * sigmai
203     else
204     myPotC = 0.0_dp
205     myDerivC = 0.0_dp
206 gezelter 1174 endif
207 chuckv 1173
208     pot_temp = epsilon * (myPot - myPotC)
209     vpair = vpair + pot_temp
210     dudr = sw * epsilon * (myDeriv - myDerivC) * sigmai
211    
212     drdx = d(1) / rij
213     drdy = d(2) / rij
214     drdz = d(3) / rij
215    
216     fx = dudr * drdx
217     fy = dudr * drdy
218     fz = dudr * drdz
219 gezelter 1174
220 chuckv 1173 #ifdef IS_MPI
221     if (do_pot) then
222     pot_Row(VDW_POT,atom1) = pot_Row(VDW_POT,atom1) + sw*pot_temp*0.5
223     pot_Col(VDW_POT,atom2) = pot_Col(VDW_POT,atom2) + sw*pot_temp*0.5
224     endif
225 gezelter 1174
226 chuckv 1173 f_Row(1,atom1) = f_Row(1,atom1) + fx
227     f_Row(2,atom1) = f_Row(2,atom1) + fy
228     f_Row(3,atom1) = f_Row(3,atom1) + fz
229 gezelter 1174
230 chuckv 1173 f_Col(1,atom2) = f_Col(1,atom2) - fx
231     f_Col(2,atom2) = f_Col(2,atom2) - fy
232     f_Col(3,atom2) = f_Col(3,atom2) - fz
233 gezelter 1174
234 chuckv 1173 #else
235     if (do_pot) pot = pot + sw*pot_temp
236 gezelter 1174
237 chuckv 1173 f(1,atom1) = f(1,atom1) + fx
238     f(2,atom1) = f(2,atom1) + fy
239     f(3,atom1) = f(3,atom1) + fz
240 gezelter 1174
241 chuckv 1173 f(1,atom2) = f(1,atom2) - fx
242     f(2,atom2) = f(2,atom2) - fy
243     f(3,atom2) = f(3,atom2) - fz
244     #endif
245    
246     #ifdef IS_MPI
247     id1 = AtomRowToGlobal(atom1)
248     id2 = AtomColToGlobal(atom2)
249     #else
250     id1 = atom1
251     id2 = atom2
252     #endif
253 gezelter 1174
254 chuckv 1173 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
255 gezelter 1174
256 chuckv 1173 fpair(1) = fpair(1) + fx
257     fpair(2) = fpair(2) + fy
258     fpair(3) = fpair(3) + fz
259 gezelter 1174
260 chuckv 1173 endif
261 gezelter 1174 return
262 chuckv 1173 end subroutine calc_mnm_lennardjones
263    
264     subroutine calc_mnm_morse(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
265     Pot, f, Do_pot, interaction_id, interaction_type)
266 gezelter 1174 integer, intent(in) :: atom1, atom2
267     real( kind = dp ), intent(in) :: rij, r2, rcut
268     real( kind = dp ) :: pot, sw, vpair
269 chuckv 1173 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
270 gezelter 1174 real( kind = dp ), intent(in), dimension(3) :: d
271 chuckv 1173 real( kind = dp ), intent(inout), dimension(3) :: fpair
272 gezelter 1174 logical, intent(in) :: do_pot
273 chuckv 1173 integer, intent(in) :: interaction_id, interaction_type
274 gezelter 1174 logical :: shiftedPot, shiftedFrc
275 chuckv 1173
276     ! Local Variables
277     real(kind=dp) :: Beta0
278     real(kind=dp) :: R0
279     real(kind=dp) :: D0
280     real(kind=dp) :: expt
281     real(kind=dp) :: expt2
282     real(kind=dp) :: expfnc
283     real(kind=dp) :: expfnc2
284     real(kind=dp) :: D_expt
285     real(kind=dp) :: D_expt2
286     real(kind=dp) :: rcos
287     real(kind=dp) :: exptC
288     real(kind=dp) :: expt2C
289     real(kind=dp) :: expfncC
290     real(kind=dp) :: expfnc2C
291     real(kind=dp) :: D_expt2C
292     real(kind=dp) :: D_exptC
293    
294     real(kind=dp) :: myPot
295     real(kind=dp) :: myPotC
296     real(kind=dp) :: myDeriv
297     real(kind=dp) :: myDerivC
298     real(kind=dp) :: pot_temp
299     real(kind=dp) :: fx,fy,fz
300     real(kind=dp) :: drdx,drdy,drdz
301     real(kind=dp) :: dudr
302     integer :: id1,id2
303 gezelter 1174
304 chuckv 1173
305 gezelter 1174 D0 = MnM_Map%interactions(interaction_id)%D0
306     R0 = MnM_Map%interactions(interaction_id)%r0
307 chuckv 1173 Beta0 = MnM_Map%interactions(interaction_id)%Beta0
308 gezelter 1174 shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot
309     shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc
310    
311     ! V(r) = D_e exp(-a(r-re)(exp(-a(r-re))-2)
312    
313     expt = -beta0*(rij - R0)
314 chuckv 1173 expfnc = exp(expt)
315 gezelter 1174 expfnc2 = expfnc*expfnc
316 chuckv 1173
317 gezelter 1174 if (shiftedPot .or. shiftedFrc) then
318     exptC = -beta0*(rcut - R0)
319     expfncC = exp(exptC)
320     expfnc2C = expfncC*expfncC
321     endif
322    
323 chuckv 1173 select case (interaction_type)
324 gezelter 1174 case (MNM_SHIFTEDMORSE)
325 chuckv 1173
326 gezelter 1174 myPot = D0 * (expfnc2 - 2.0_dp * expfnc)
327     myDeriv = 2.0*D0*beta0*(expfnc - expfnc2)
328 chuckv 1173
329 gezelter 1174 if (shiftedPot) then
330     myPotC = D0 * (expfnc2C - 2.0_dp * expfncC)
331     myDerivC = 0.0_dp
332     elseif (shiftedFrc) then
333     myPotC = D0 * (expfnc2C - 2.0_dp * expfncC)
334     myDerivC = 2.0*D0*beta0*(expfnc2C - expfnc2C)
335     myPotC = myPotC + myDerivC * (rij - rcut)
336     else
337     myPotC = 0.0_dp
338     myDerivC = 0.0_dp
339     endif
340 chuckv 1173
341     case (MNM_REPULSIVEMORSE)
342    
343 gezelter 1174 myPot = D0 * expfnc2
344     myDeriv = -2.0_dp * D0 * beta0 * expfnc2
345 chuckv 1173
346 gezelter 1174 if (shiftedPot) then
347     myPotC = D0 * expfnc2C
348     myDerivC = 0.0_dp
349     elseif (shiftedFrc) then
350     myPotC = D0 * expfnc2C
351     myDerivC = -2.0_dp * D0* beta0 * expfnc2C
352     myPotC = myPotC + myDerivC * (rij - rcut)
353     else
354     myPotC = 0.0_dp
355     myDerivC = 0.0_dp
356     endif
357 chuckv 1173 end select
358    
359     pot_temp = (myPot - myPotC)
360     vpair = vpair + pot_temp
361     dudr = sw * (myDeriv - myDerivC)
362    
363     drdx = d(1) / rij
364     drdy = d(2) / rij
365     drdz = d(3) / rij
366    
367     fx = dudr * drdx
368     fy = dudr * drdy
369     fz = dudr * drdz
370    
371     #ifdef IS_MPI
372     if (do_pot) then
373     pot_Row(VDW_POT,atom1) = pot_Row(VDW_POT,atom1) + sw*pot_temp*0.5
374     pot_Col(VDW_POT,atom2) = pot_Col(VDW_POT,atom2) + sw*pot_temp*0.5
375     endif
376    
377     f_Row(1,atom1) = f_Row(1,atom1) + fx
378     f_Row(2,atom1) = f_Row(2,atom1) + fy
379     f_Row(3,atom1) = f_Row(3,atom1) + fz
380    
381     f_Col(1,atom2) = f_Col(1,atom2) - fx
382     f_Col(2,atom2) = f_Col(2,atom2) - fy
383     f_Col(3,atom2) = f_Col(3,atom2) - fz
384    
385     #else
386     if (do_pot) pot = pot + sw*pot_temp
387    
388     f(1,atom1) = f(1,atom1) + fx
389     f(2,atom1) = f(2,atom1) + fy
390     f(3,atom1) = f(3,atom1) + fz
391    
392     f(1,atom2) = f(1,atom2) - fx
393     f(2,atom2) = f(2,atom2) - fy
394     f(3,atom2) = f(3,atom2) - fz
395     #endif
396    
397     #ifdef IS_MPI
398     id1 = AtomRowToGlobal(atom1)
399     id2 = AtomColToGlobal(atom2)
400     #else
401     id1 = atom1
402     id2 = atom2
403     #endif
404    
405     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
406    
407     fpair(1) = fpair(1) + fx
408     fpair(2) = fpair(2) + fy
409     fpair(3) = fpair(3) + fz
410    
411     endif
412    
413     return
414     end subroutine calc_mnm_morse
415    
416     subroutine calc_mnm_maw(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
417 gezelter 1174 Pot, A, F,t, Do_pot, interaction_id)
418     integer, intent(in) :: atom1, atom2
419     real( kind = dp ), intent(in) :: rij, r2, rcut
420     real( kind = dp ) :: pot, sw, vpair
421 chuckv 1173 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
422 gezelter 1174 real (kind=dp),intent(in), dimension(9,nLocal) :: A
423 chuckv 1173 real (kind=dp),intent(inout), dimension(3,nLocal) :: t
424    
425 gezelter 1174 real( kind = dp ), intent(in), dimension(3) :: d
426 chuckv 1173 real( kind = dp ), intent(inout), dimension(3) :: fpair
427 gezelter 1174 logical, intent(in) :: do_pot
428 chuckv 1173
429     integer, intent(in) :: interaction_id
430    
431 gezelter 1238 real(kind = dp) :: D0, R0, beta0, ca1, cb1, pot_temp
432 gezelter 1174 real(kind = dp) :: expt0, expfnc0, expfnc02
433     real(kind = dp) :: exptH, expfncH, expfncH2
434 chuckv 1231 real(kind = dp) :: x, y, z, x2, y2, z2, r3, r4
435     real(kind = dp) :: drdx, drdy, drdz
436     real(kind = dp) :: dvdx, dvdy, dvdz
437 gezelter 1238 real(kind = dp) :: Vang, dVangdx, dVangdy, dVangdz
438     real(kind = dp) :: dVangdux, dVangduy, dVangduz
439 chuckv 1231 real(kind = dp) :: dVmorsedr
440 chuckv 1229 real(kind = dp) :: Vmorse, dVmorsedx, dVmorsedy, dVmorsedz
441 gezelter 1238 real(kind = dp) :: a1, b1, s
442     real(kind = dp) :: da1dx, da1dy, da1dz, da1dux, da1duy, da1duz
443     real(kind = dp) :: db1dx, db1dy, db1dz, db1dux, db1duy, db1duz
444 gezelter 1174 real(kind = dp) :: fx, fy, fz, tx, ty, tz, fxl, fyl, fzl
445 gezelter 1238 real(kind = dp), parameter :: st = sqrt(3.0_dp)
446 gezelter 1174 integer :: atid1, atid2, id1, id2
447     logical :: shiftedPot, shiftedFrc
448    
449     #ifdef IS_MPI
450     atid1 = atid_Row(atom1)
451     atid2 = atid_Col(atom2)
452    
453     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
454     ! rotate the inter-particle separation into the two different
455     ! body-fixed coordinate systems:
456    
457     x = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
458     y = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
459     z = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
460     else
461     ! negative sign because this is the vector from j to i:
462    
463     x = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
464     y = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
465     z = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
466     endif
467    
468     #else
469     atid1 = atid(atom1)
470     atid2 = atid(atom2)
471 chuckv 1175
472 gezelter 1174 if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
473     ! rotate the inter-particle separation into the two different
474     ! body-fixed coordinate systems:
475    
476     x = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
477     y = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
478     z = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
479     else
480     ! negative sign because this is the vector from j to i:
481    
482     x = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
483     y = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
484     z = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
485     endif
486    
487     #endif
488 chuckv 1231
489 gezelter 1174 D0 = MnM_Map%interactions(interaction_id)%D0
490     R0 = MnM_Map%interactions(interaction_id)%r0
491 chuckv 1229 beta0 = MnM_Map%interactions(interaction_id)%beta0
492 gezelter 1238 ca1 = MnM_Map%interactions(interaction_id)%ca1
493     cb1 = MnM_Map%interactions(interaction_id)%cb1
494 chuckv 1175
495 gezelter 1174 shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot
496     shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc
497    
498     expt0 = -beta0*(rij - R0)
499     expfnc0 = exp(expt0)
500     expfnc02 = expfnc0*expfnc0
501 chuckv 1175
502 gezelter 1174 !!$ if (shiftedPot .or. shiftedFrc) then
503     !!$ exptC0 = -beta0*(rcut - R0)
504     !!$ expfncC0 = exp(exptC0)
505     !!$ expfncC02 = expfncC0*expfncC0
506     !!$ exptCH = -betaH*(rcut - R0)
507     !!$ expfncCH = exp(exptCH)
508     !!$ expfncCH2 = expfncCH*expfncCH
509     !!$ endif
510    
511 chuckv 1231 drdx = x / rij
512     drdy = y / rij
513     drdz = z / rij
514    
515 gezelter 1174 x2 = x*x
516     y2 = y*y
517     z2 = z*z
518 chuckv 1175 r3 = r2*rij
519 chuckv 1231 r4 = r2*r2
520 chuckv 1176
521 chuckv 1231 Vmorse = D0 * (expfnc02 - 2.0_dp * expfnc0)
522 chuckv 1176
523 gezelter 1238 ! angular modulation of morse part of potential to approximate
524     ! the squares of the two HOMO lone pair orbitals in water:
525 chuckv 1245 !********************** old form*************************
526 gezelter 1238 ! s = 1 / (4 pi)
527 chuckv 1242 ! a1 = (s - pz)^2 = (1 - sqrt(3)*cos(theta))^2 / (4 pi)
528 gezelter 1238 ! b1 = px^2 = 3 * (sin(theta)*cos(phi))^2 / (4 pi)
529 chuckv 1245 !********************** old form*************************
530 gezelter 1238 ! we'll leave out the 4 pi for now
531 chuckv 1245
532     ! new functional form just using the p orbitals.
533     ! Vmorse(r)*[a*p_x + b p_z + (1-a-b)]
534     ! which is
535     ! Vmorse(r)*[a sin^2(theta) cos^2(phi) + b cos(theta) + (1-a-b)]
536     ! Vmorse(r)*[a*x2/r2 + b*z/r + (1-a-b)]
537 gezelter 1238
538 chuckv 1245
539    
540 gezelter 1238 s = 1.0_dp
541 chuckv 1245 ! a1 = (1.0_dp - st * z / rij )**2
542     ! b1 = 3.0_dp * x2 / r2
543 gezelter 1238
544 chuckv 1245 ! Vang = s + ca1 * a1 + cb1 * b1
545    
546     Vang = ca1 * x2/r2 + cb1 * z/rij + (0.8_dp-ca1/3.0_dp)
547    
548 chuckv 1229 pot_temp = Vmorse*Vang
549 chuckv 1228
550 gezelter 1174 vpair = vpair + pot_temp
551 chuckv 1175
552 gezelter 1174 if (do_pot) then
553     #ifdef IS_MPI
554     pot_row(VDW_POT,atom1) = pot_row(VDW_POT,atom1) + 0.5_dp*pot_temp*sw
555     pot_col(VDW_POT,atom2) = pot_col(VDW_POT,atom2) + 0.5_dp*pot_temp*sw
556     #else
557     pot = pot + pot_temp*sw
558     #endif
559     endif
560    
561 chuckv 1229 dVmorsedr = 2.0_dp*D0*beta0*(expfnc0 - expfnc02)
562 chuckv 1231
563 chuckv 1229 dVmorsedx = dVmorsedr * drdx
564     dVmorsedy = dVmorsedr * drdy
565     dVmorsedz = dVmorsedr * drdz
566 gezelter 1238
567 chuckv 1245 !da1dx = 2.0_dp * st * x * z / r3 - 6.0_dp * x * z2 / r4
568     !da1dy = 2.0_dp * st * y * z / r3 - 6.0_dp * y * z2 / r4
569     !da1dz = 2.0_dp * st * (x2 + y2) * (st * z - rij ) / r4
570 chuckv 1231
571 chuckv 1245 !db1dx = 6.0_dp * x * (1.0_dp - x2 / r2) / r2
572     !db1dy = -6.0_dp * x2 * y / r4
573     !db1dz = -6.0_dp * x2 * z / r4
574 gezelter 1238
575 chuckv 1245 !dVangdx = ca1 * da1dx + cb1 * db1dx
576     !dVangdy = ca1 * da1dy + cb1 * db1dy
577     !dVangdz = ca1 * da1dz + cb1 * db1dz
578 chuckv 1228
579 chuckv 1245 dVangdx = -2.0*ca1*x2*x/r4 + 2.0*ca1*x/r2 - cb1*x*z/r3
580     dVangdy = -2.0*ca1*x2*y/r4 - cb1*z*y/r3
581     dVangdz = -2.0*ca1*x2*z/r4 + cb1/rij - cb1*z2 /r3
582    
583 chuckv 1231 ! chain rule to put these back on x, y, z
584 chuckv 1229 dvdx = Vang * dVmorsedx + Vmorse * dVangdx
585     dvdy = Vang * dVmorsedy + Vmorse * dVangdy
586     dvdz = Vang * dVmorsedz + Vmorse * dVangdz
587 gezelter 1174
588 chuckv 1231 ! Torques for Vang using method of Price:
589     ! S. L. Price, A. J. Stone, and M. Alderton, Mol. Phys. 52, 987 (1984).
590 gezelter 1238
591 chuckv 1245 !da1dux = 6.0_dp * y * z / r2 - 2.0_dp * st * y / rij
592     !da1duy = -6.0_dp * x * z / r2 + 2.0_dp * st * x / rij
593     !da1duz = 0.0_dp
594 gezelter 1238
595 chuckv 1245 !db1dux = 0.0_dp
596     !db1duy = 6.0_dp * x * z / r2
597     !db1duz = -6.0_dp * y * x / r2
598 gezelter 1238
599 chuckv 1245 !dVangdux = ca1 * da1dux + cb1 * db1dux
600     !dVangduy = ca1 * da1duy + cb1 * db1duy
601     !dVangduz = ca1 * da1duz + cb1 * db1duz
602 chuckv 1231
603 chuckv 1245 dVangdux = cb1*y/rij
604     dVangduy = 2.0*ca1*x*z/r2 - cb1*x/rij
605     dVangduz = -2.0*ca1*y*x/r2
606    
607 chuckv 1231 ! do the torques first since they are easy:
608     ! remember that these are still in the body fixed axes
609    
610     tx = Vmorse * dVangdux * sw
611     ty = Vmorse * dVangduy * sw
612     tz = Vmorse * dVangduz * sw
613 gezelter 1174
614     ! go back to lab frame using transpose of rotation matrix:
615    
616     #ifdef IS_MPI
617     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
618     t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*tx + &
619     a_Row(4,atom1)*ty + a_Row(7,atom1)*tz
620     t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*tx + &
621     a_Row(5,atom1)*ty + a_Row(8,atom1)*tz
622     t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*tx + &
623     a_Row(6,atom1)*ty + a_Row(9,atom1)*tz
624     else
625     t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*tx + &
626     a_Col(4,atom2)*ty + a_Col(7,atom2)*tz
627     t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*tx + &
628     a_Col(5,atom2)*ty + a_Col(8,atom2)*tz
629     t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*tx + &
630     a_Col(6,atom2)*ty + a_Col(9,atom2)*tz
631     endif
632     #else
633     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
634     t(1,atom1) = t(1,atom1) + a(1,atom1)*tx + a(4,atom1)*ty + &
635     a(7,atom1)*tz
636     t(2,atom1) = t(2,atom1) + a(2,atom1)*tx + a(5,atom1)*ty + &
637     a(8,atom1)*tz
638     t(3,atom1) = t(3,atom1) + a(3,atom1)*tx + a(6,atom1)*ty + &
639     a(9,atom1)*tz
640     else
641     t(1,atom2) = t(1,atom2) + a(1,atom2)*tx + a(4,atom2)*ty + &
642     a(7,atom2)*tz
643     t(2,atom2) = t(2,atom2) + a(2,atom2)*tx + a(5,atom2)*ty + &
644     a(8,atom2)*tz
645     t(3,atom2) = t(3,atom2) + a(3,atom2)*tx + a(6,atom2)*ty + &
646     a(9,atom2)*tz
647     endif
648     #endif
649 chuckv 1231 ! Now, on to the forces (still in body frame of water)
650 gezelter 1174
651 chuckv 1229 fx = dvdx * sw
652     fy = dvdy * sw
653     fz = dvdz * sw
654 gezelter 1174
655     ! rotate the terms back into the lab frame:
656    
657     #ifdef IS_MPI
658     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
659     fxl = a_Row(1,atom1)*fx + a_Row(4,atom1)*fy + a_Row(7,atom1)*fz
660     fyl = a_Row(2,atom1)*fx + a_Row(5,atom1)*fy + a_Row(8,atom1)*fz
661     fzl = a_Row(3,atom1)*fx + a_Row(6,atom1)*fy + a_Row(9,atom1)*fz
662     else
663 chuckv 1234 ! negative sign because this is the vector from j to i:
664     fxl = -(a_Col(1,atom2)*fx + a_Col(4,atom2)*fy + a_Col(7,atom2)*fz)
665     fyl = -(a_Col(2,atom2)*fx + a_Col(5,atom2)*fy + a_Col(8,atom2)*fz)
666     fzl = -(a_Col(3,atom2)*fx + a_Col(6,atom2)*fy + a_Col(9,atom2)*fz)
667 gezelter 1174 endif
668     f_Row(1,atom1) = f_Row(1,atom1) + fxl
669     f_Row(2,atom1) = f_Row(2,atom1) + fyl
670     f_Row(3,atom1) = f_Row(3,atom1) + fzl
671    
672     f_Col(1,atom2) = f_Col(1,atom2) - fxl
673     f_Col(2,atom2) = f_Col(2,atom2) - fyl
674     f_Col(3,atom2) = f_Col(3,atom2) - fzl
675     #else
676     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
677     fxl = a(1,atom1)*fx + a(4,atom1)*fy + a(7,atom1)*fz
678     fyl = a(2,atom1)*fx + a(5,atom1)*fy + a(8,atom1)*fz
679     fzl = a(3,atom1)*fx + a(6,atom1)*fy + a(9,atom1)*fz
680     else
681 chuckv 1234 ! negative sign because this is the vector from j to i:
682     fxl = -(a(1,atom2)*fx + a(4,atom2)*fy + a(7,atom2)*fz)
683     fyl = -(a(2,atom2)*fx + a(5,atom2)*fy + a(8,atom2)*fz)
684     fzl = -(a(3,atom2)*fx + a(6,atom2)*fy + a(9,atom2)*fz)
685 gezelter 1174 endif
686     f(1,atom1) = f(1,atom1) + fxl
687     f(2,atom1) = f(2,atom1) + fyl
688     f(3,atom1) = f(3,atom1) + fzl
689    
690     f(1,atom2) = f(1,atom2) - fxl
691     f(2,atom2) = f(2,atom2) - fyl
692     f(3,atom2) = f(3,atom2) - fzl
693     #endif
694    
695     #ifdef IS_MPI
696     id1 = AtomRowToGlobal(atom1)
697     id2 = AtomColToGlobal(atom2)
698     #else
699     id1 = atom1
700     id2 = atom2
701     #endif
702    
703     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
704    
705     fpair(1) = fpair(1) + fxl
706     fpair(2) = fpair(2) + fyl
707     fpair(3) = fpair(3) + fzl
708    
709     endif
710    
711     return
712 chuckv 1173 end subroutine calc_mnm_maw
713    
714    
715 chuckv 1168 subroutine setMnMDefaultCutoff(thisRcut, shiftedPot, shiftedFrc)
716     real(kind=dp), intent(in) :: thisRcut
717     logical, intent(in) :: shiftedPot
718     logical, intent(in) :: shiftedFrc
719     integer i, nInteractions
720     defaultCutoff = thisRcut
721     defaultShiftPot = shiftedPot
722     defaultShiftFrc = shiftedFrc
723    
724 xsun 1178 if (associated(MnM_Map)) then
725     if(MnM_Map%interactionCount /= 0) then
726     nInteractions = MnM_Map%interactionCount
727 chuckv 1168
728 xsun 1178 do i = 1, nInteractions
729     MnM_Map%interactions(i)%shiftedPot = shiftedPot
730     MnM_Map%interactions(i)%shiftedFrc = shiftedFrc
731     MnM_Map%interactions(i)%rCut = thisRcut
732     MnM_Map%interactions(i)%rCutWasSet = .true.
733     enddo
734     end if
735     end if
736 chuckv 1168
737     end subroutine setMnMDefaultCutoff
738    
739     subroutine copyAllData(v1, v2)
740     type(MnMinteractionMap), pointer :: v1
741     type(MnMinteractionMap), pointer :: v2
742     integer :: i, j
743    
744     do i = 1, v1%interactionCount
745     v2%interactions(i) = v1%interactions(i)
746     enddo
747    
748     v2%interactionCount = v1%interactionCount
749     return
750     end subroutine copyAllData
751    
752     subroutine addInteraction(myInteraction)
753     type(MNMtype) :: myInteraction
754     type(MnMinteraction) :: nt
755     integer :: id
756 chuckv 1175
757 chuckv 1168 nt%interaction_type = myInteraction%MNMInteractionType
758 chuckv 1175 nt%metal_atid = &
759     getFirstMatchingElement(atypes, "c_ident", myInteraction%metal_atid)
760     nt%nonmetal_atid = &
761     getFirstMatchingElement(atypes, "c_ident", myInteraction%nonmetal_atid)
762 gezelter 1238
763 chuckv 1168 select case (nt%interaction_type)
764     case (MNM_LENNARDJONES)
765     nt%sigma = myInteraction%sigma
766     nt%epsilon = myInteraction%epsilon
767     case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE)
768     nt%R0 = myInteraction%R0
769     nt%D0 = myInteraction%D0
770     nt%beta0 = myInteraction%beta0
771     case(MNM_MAW)
772     nt%R0 = myInteraction%R0
773     nt%D0 = myInteraction%D0
774     nt%beta0 = myInteraction%beta0
775 gezelter 1238 nt%ca1 = myInteraction%ca1
776     nt%cb1 = myInteraction%cb1
777 chuckv 1168 case default
778 chuckv 1173 call handleError("MNM", "Unknown Interaction type")
779 chuckv 1168 end select
780    
781     if (.not. associated(MnM_Map)) then
782     call ensureCapacityHelper(MnM_Map, 1)
783     else
784     call ensureCapacityHelper(MnM_Map, MnM_Map%interactionCount + 1)
785     end if
786    
787     MnM_Map%interactionCount = MnM_Map%interactionCount + 1
788     id = MnM_Map%interactionCount
789     MnM_Map%interactions(id) = nt
790     end subroutine addInteraction
791    
792     subroutine ensureCapacityHelper(this, minCapacity)
793     type(MnMinteractionMap), pointer :: this, that
794     integer, intent(in) :: minCapacity
795     integer :: oldCapacity
796     integer :: newCapacity
797     logical :: resizeFlag
798    
799     resizeFlag = .false.
800    
801     ! first time: allocate a new vector with default size
802    
803     if (.not. associated(this)) then
804     this => MnMinitialize(minCapacity, 0)
805     endif
806    
807     oldCapacity = size(this%interactions)
808    
809     if (minCapacity > oldCapacity) then
810     if (this%capacityIncrement .gt. 0) then
811     newCapacity = oldCapacity + this%capacityIncrement
812     else
813     newCapacity = oldCapacity * 2
814     endif
815     if (newCapacity .lt. minCapacity) then
816     newCapacity = minCapacity
817     endif
818     resizeFlag = .true.
819     else
820     newCapacity = oldCapacity
821     endif
822    
823     if (resizeFlag) then
824     that => MnMinitialize(newCapacity, this%capacityIncrement)
825     call copyAllData(this, that)
826     this => MnMdestroy(this)
827     this => that
828     endif
829     end subroutine ensureCapacityHelper
830    
831     function MnMinitialize(cap, capinc) result(this)
832     integer, intent(in) :: cap, capinc
833     integer :: error
834     type(MnMinteractionMap), pointer :: this
835    
836     nullify(this)
837    
838     if (cap < 0) then
839     write(*,*) 'Bogus Capacity:', cap
840     return
841     endif
842     allocate(this,stat=error)
843     if ( error /= 0 ) then
844     write(*,*) 'Could not allocate MnMinteractionMap!'
845     return
846     end if
847    
848     this%initialCapacity = cap
849     this%capacityIncrement = capinc
850    
851     allocate(this%interactions(this%initialCapacity), stat=error)
852     if(error /= 0) write(*,*) 'Could not allocate MnMinteraction!'
853    
854     end function MnMinitialize
855    
856 chuckv 1175 subroutine createInteractionLookup(this)
857     type (MnMInteractionMap),pointer :: this
858 chuckv 1168 integer :: biggestAtid, i, metal_atid, nonmetal_atid, error
859    
860     biggestAtid=-1
861     do i = 1, this%interactionCount
862     metal_atid = this%interactions(i)%metal_atid
863     nonmetal_atid = this%interactions(i)%nonmetal_atid
864    
865     if (metal_atid .gt. biggestAtid) biggestAtid = metal_atid
866     if (nonmetal_atid .gt. biggestAtid) biggestAtid = nonmetal_atid
867     enddo
868 chuckv 1175
869 chuckv 1168 allocate(MnMinteractionLookup(biggestAtid,biggestAtid), stat=error)
870     if (error /= 0) write(*,*) 'Could not allocate MnMinteractionLookup'
871    
872     do i = 1, this%interactionCount
873     metal_atid = this%interactions(i)%metal_atid
874     nonmetal_atid = this%interactions(i)%nonmetal_atid
875    
876     MnMinteractionLookup(metal_atid, nonmetal_atid) = i
877     MnMinteractionLookup(nonmetal_atid, metal_atid) = i
878     enddo
879     end subroutine createInteractionLookup
880    
881     function MnMdestroy(this) result(null_this)
882     logical :: done
883     type(MnMinteractionMap), pointer :: this
884     type(MnMinteractionMap), pointer :: null_this
885    
886     if (.not. associated(this)) then
887     null_this => null()
888     return
889     end if
890    
891     !! Walk down the list and deallocate each of the map's components
892     if(associated(this%interactions)) then
893     deallocate(this%interactions)
894     this%interactions=>null()
895     endif
896     deallocate(this)
897     this => null()
898     null_this => null()
899     end function MnMdestroy
900    
901     subroutine deleteInteractions()
902     MnM_Map => MnMdestroy(MnM_Map)
903     return
904     end subroutine deleteInteractions
905    
906 chuckv 1173 subroutine getLJfunc(r, myPot, myDeriv)
907    
908     real(kind=dp), intent(in) :: r
909     real(kind=dp), intent(inout) :: myPot, myDeriv
910     real(kind=dp) :: ri, ri2, ri6, ri7, ri12, ri13
911     real(kind=dp) :: a, b, c, d, dx
912     integer :: j
913    
914     ri = 1.0_DP / r
915     ri2 = ri*ri
916     ri6 = ri2*ri2*ri2
917     ri7 = ri6*ri
918     ri12 = ri6*ri6
919     ri13 = ri12*ri
920    
921     myPot = 4.0_DP * (ri12 - ri6)
922     myDeriv = 24.0_DP * (ri7 - 2.0_DP * ri13)
923    
924     return
925     end subroutine getLJfunc
926    
927     subroutine getSoftFunc(r, myPot, myDeriv)
928    
929     real(kind=dp), intent(in) :: r
930     real(kind=dp), intent(inout) :: myPot, myDeriv
931     real(kind=dp) :: ri, ri2, ri6, ri7
932     real(kind=dp) :: a, b, c, d, dx
933     integer :: j
934    
935     ri = 1.0_DP / r
936     ri2 = ri*ri
937     ri6 = ri2*ri2*ri2
938     ri7 = ri6*ri
939     myPot = 4.0_DP * (ri6)
940     myDeriv = - 24.0_DP * ri7
941    
942     return
943     end subroutine getSoftFunc
944 chuckv 1168 end module MetalNonMetal