ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 1228
Committed: Thu Feb 14 21:37:05 2008 UTC (17 years, 5 months ago) by chuckv
Original Path: trunk/src/UseTheForce/DarkSide/MetalNonMetal.F90
File size: 30934 byte(s)
Log Message:
Changes to simparalell to remove MPI stuff.

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 1228 !! @version $Id: MetalNonMetal.F90,v 1.7 2008-02-14 21:37:05 chuckv Exp $, $Date: 2008-02-14 21:37:05 $, $Name: not supported by cvs2svn $, $Revision: 1.7 $
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     real(kind=DP), parameter :: sqrt3_2 = sqrt(3.0_dp)/2.0_dp
72     real(kind=DP), parameter :: threefourth = 3.0_dp/4.0_dp
73     real(kind=DP), parameter :: sqrt3_4 = sqrt(3.0_dp)/4.0_dp
74    
75 chuckv 1168 logical, save :: defaultShiftPot = .false.
76     logical, save :: defaultShiftFrc = .false.
77     logical, save :: haveDefaultCutoff = .false.
78    
79     type :: MnMinteraction
80     integer :: metal_atid
81     integer :: nonmetal_atid
82     integer :: interaction_type
83     real(kind=dp) :: R0
84     real(kind=dp) :: D0
85     real(kind=dp) :: beta0
86     real(kind=dp) :: betaH
87     real(kind=dp) :: alpha
88     real(kind=dp) :: gamma
89     real(kind=dp) :: sigma
90     real(kind=dp) :: epsilon
91     real(kind=dp) :: rCut = 0.0_dp
92     logical :: rCutWasSet = .false.
93     logical :: shiftedPot = .false.
94     logical :: shiftedFrc = .false.
95     end type MnMinteraction
96    
97     type :: MnMinteractionMap
98     PRIVATE
99     integer :: initialCapacity = 10
100     integer :: capacityIncrement = 0
101     integer :: interactionCount = 0
102     type(MnMinteraction), pointer :: interactions(:) => null()
103     end type MnMinteractionMap
104    
105 xsun 1178 type (MnMInteractionMap), pointer :: MnM_Map
106 chuckv 1168
107     integer, allocatable, dimension(:,:) :: MnMinteractionLookup
108    
109     public :: setMnMDefaultCutoff
110     public :: addInteraction
111     public :: deleteInteractions
112     public :: MNMtype
113 chuckv 1173 public :: do_mnm_pair
114 chuckv 1168
115     contains
116    
117    
118 chuckv 1173 subroutine do_mnm_pair(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
119     Pot, A, F,t, Do_pot)
120 gezelter 1174 integer, intent(in) :: atom1, atom2
121 chuckv 1173 integer :: atid1, atid2, ljt1, ljt2
122 gezelter 1174 real( kind = dp ), intent(in) :: rij, r2, rcut
123     real( kind = dp ) :: pot, sw, vpair
124 chuckv 1173 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
125 gezelter 1174 real (kind=dp), intent(in), dimension(9,nLocal) :: A
126     real (kind=dp), intent(inout), dimension(3,nLocal) :: t
127     real( kind = dp ), intent(in), dimension(3) :: d
128 chuckv 1173 real( kind = dp ), intent(inout), dimension(3) :: fpair
129 gezelter 1174 logical, intent(in) :: do_pot
130 chuckv 1168
131 chuckv 1173 integer :: interaction_id
132     integer :: interaction_type
133    
134     #ifdef IS_MPI
135     atid1 = atid_Row(atom1)
136     atid2 = atid_Col(atom2)
137     #else
138     atid1 = atid(atom1)
139     atid2 = atid(atom2)
140     #endif
141    
142 chuckv 1175 if(.not.haveInteractionLookup) then
143     call createInteractionLookup(MnM_MAP)
144     haveInteractionLookup =.true.
145     end if
146    
147 chuckv 1173 interaction_id = MnMinteractionLookup(atid1, atid2)
148     interaction_type = MnM_Map%interactions(interaction_id)%interaction_type
149    
150     select case (interaction_type)
151     case (MNM_LENNARDJONES)
152 gezelter 1174 call calc_mnm_lennardjones(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, &
153     Fpair, Pot, F, Do_pot, interaction_id)
154 chuckv 1173 case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE)
155     call calc_mnm_morse(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
156 gezelter 1174 Pot, F, Do_pot, interaction_id, interaction_type)
157 chuckv 1173 case(MNM_MAW)
158     call calc_mnm_maw(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
159 gezelter 1174 Pot,A, F,t, Do_pot, interaction_id)
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 1174 real(kind = dp) :: D0, R0, beta0, betaH, alpha, gamma, pot_temp
432     real(kind = dp) :: expt0, expfnc0, expfnc02
433     real(kind = dp) :: exptH, expfncH, expfncH2
434 chuckv 1228 real(kind = dp) :: x, y, z, x2, y2, z2, r3, proj, proj3, st2,st,s2t
435 gezelter 1174 real(kind = dp) :: drdx, drdy, drdz, drdux, drduy, drduz
436 chuckv 1228 real(kind = dp) :: ct,ct2,c2t, dctdx, dctdy, dctdz, dctdux, dctduy, dctduz
437     real(kind = dp) :: sp,c2p,sp2,dspdx, dspdy, dspdz, dspdux, dspduy, dspduz
438 gezelter 1174 real(kind = dp) :: dvdx, dvdy, dvdz, dvdux, dvduy, dvduz
439 chuckv 1228 real(kind = dp) :: maw, dmawdr, dmawdct, dmawdsp, dmawdt, dmawdst
440     real(kind = dp) :: hpart, vmorse
441 gezelter 1174 real(kind = dp) :: fx, fy, fz, tx, ty, tz, fxl, fyl, fzl
442     integer :: atid1, atid2, id1, id2
443     logical :: shiftedPot, shiftedFrc
444    
445 chuckv 1175
446    
447    
448 gezelter 1174 #ifdef IS_MPI
449     atid1 = atid_Row(atom1)
450     atid2 = atid_Col(atom2)
451    
452     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
453     ! rotate the inter-particle separation into the two different
454     ! body-fixed coordinate systems:
455    
456     x = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
457     y = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
458     z = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
459     else
460     ! negative sign because this is the vector from j to i:
461    
462     x = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
463     y = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
464     z = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
465     endif
466    
467     #else
468     atid1 = atid(atom1)
469     atid2 = atid(atom2)
470 chuckv 1175
471 gezelter 1174 if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
472     ! rotate the inter-particle separation into the two different
473     ! body-fixed coordinate systems:
474    
475     x = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
476     y = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
477     z = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
478     else
479     ! negative sign because this is the vector from j to i:
480    
481     x = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
482     y = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
483     z = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
484     endif
485    
486     #endif
487    
488 chuckv 1175
489 gezelter 1174 D0 = MnM_Map%interactions(interaction_id)%D0
490     R0 = MnM_Map%interactions(interaction_id)%r0
491     beta0 = MnM_Map%interactions(interaction_id)%beta0
492     betaH = MnM_Map%interactions(interaction_id)%betaH
493     alpha = MnM_Map%interactions(interaction_id)%alpha
494     gamma = MnM_Map%interactions(interaction_id)%gamma
495 chuckv 1175
496 gezelter 1174
497     shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot
498     shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc
499    
500     expt0 = -beta0*(rij - R0)
501     expfnc0 = exp(expt0)
502     expfnc02 = expfnc0*expfnc0
503    
504     exptH = -betaH*(rij - R0)
505     expfncH = exp(exptH)
506     expfncH2 = expfncH*expfncH
507    
508 chuckv 1175
509    
510 gezelter 1174 !!$ if (shiftedPot .or. shiftedFrc) then
511     !!$ exptC0 = -beta0*(rcut - R0)
512     !!$ expfncC0 = exp(exptC0)
513     !!$ expfncC02 = expfncC0*expfncC0
514     !!$ exptCH = -betaH*(rcut - R0)
515     !!$ expfncCH = exp(exptCH)
516     !!$ expfncCH2 = expfncCH*expfncCH
517     !!$ endif
518    
519     drdx = -d(1) / rij
520     drdy = -d(2) / rij
521     drdz = -d(3) / rij
522     drdux = 0.0_dp
523     drduy = 0.0_dp
524     drduz = 0.0_dp
525    
526     x2 = x*x
527     y2 = y*y
528     z2 = z*z
529 chuckv 1175 r3 = r2*rij
530 gezelter 1174 ct = z / rij
531 chuckv 1228 ct2 = ct * ct
532     ! cos(2 theta)
533     c2t = 2.0_dp * ct2 - 1.0_dp
534 gezelter 1174
535     if (ct .gt. 1.0_dp) ct = 1.0_dp
536     if (ct .lt. -1.0_dp) ct = -1.0_dp
537 chuckv 1176
538     ! These derivatives can be obtained by using
539     ! cos(theta) = \hat{u} . \vec{r} / r
540     ! where \hat{u} is the body-fixed unit frame for the water molecule,
541     ! and \vec{r} is the vector to the metal atom.
542     !
543     ! the derivatives wrt \vec{r} are:
544     ! dcostheta/d\vec{r} = - costheta \vec{r} / r^2 + \hat{u} / r
545     ! the molecular frame for each water has u = {0, 0, 1}, so these:
546     !
547     ! dctdx = - x * z / r3 + ux / rij
548     ! dctdy = - y * z / r3 + uy / rij
549     ! dctdz = - z2 / r3 + uz / rij
550     !
551     ! become:
552     !
553 gezelter 1174 dctdx = - z * x / r3
554     dctdy = - z * y / r3
555     dctdz = 1.0_dp / rij - z2 / r3
556 chuckv 1176
557     dctdux = x / rij
558     dctduy = y / rij
559     dctduz = z / rij
560    
561     ! dctdux = y / rij ! - (z * x2) / r3
562     ! dctduy = -x / rij !- (z * y2) / r3
563     ! dctduz = 0.0_dp ! z / rij - (z2 * z) / r3
564 chuckv 1175
565 gezelter 1174 ! this is an attempt to try to truncate the singularity when
566     ! sin(theta) is near 0.0:
567    
568 chuckv 1228 st2 = 1.0_dp - ct2
569    
570    
571 gezelter 1174 if (abs(st2) .lt. 1.0e-12_dp) then
572     proj = sqrt(rij * 1.0e-12_dp)
573     !!$ dcpdx = 1.0_dp / proj
574     !!$ dcpdy = 0.0_dp
575     !!$ dcpdux = x / proj
576     !!$ dcpduy = 0.0_dp
577     dspdx = 0.0_dp
578     dspdy = 1.0_dp / proj
579     dspdux = 0.0_dp
580     dspduy = y / proj
581     else
582     proj = sqrt(x2 + y2)
583     proj3 = proj*proj*proj
584     !!$ dcpdx = 1.0_dp / proj - x2 / proj3
585     !!$ dcpdy = - x * y / proj3
586     !!$ dcpdux = x / proj - (x2 * x) / proj3
587     !!$ dcpduy = - (x * y2) / proj3
588     dspdx = - x * y / proj3
589     dspdy = 1.0_dp / proj - y2 / proj3
590     dspdux = - (y * x2) / proj3
591     dspduy = y / proj - (y2 * y) / proj3
592     endif
593    
594     !!$ cp = x / proj
595     !!$ dcpdz = 0.0_dp
596     !!$ dcpduz = 0.0_dp
597    
598     sp = y / proj
599     dspdz = 0.0_dp
600     dspduz = 0.0_dp
601 chuckv 1175
602 chuckv 1228 sp2 = sp * sp
603     c2p = 1.0_dp - 2.0_dp * sp2
604 gezelter 1174
605 chuckv 1228
606     ! V(r) = D_e exp(-a(r-re)(exp(-a(r-re))-2)
607     vmorse =D0 * (expfnc02 - 2.0_dp * expfnc0)
608     ! angular modulation of morse part of potential to approximate a sp3 orbital
609     ! hpart = (sqrt(3.0)*cos(theta)/2.0 - 3.0*cos(2.0*phi)*sin^2(theta)/4.0);
610    
611     hpart = (sqrt3_2*ct - threefourth*c2p*st2)
612    
613     pot_temp = vmorse* (1.0_dp - alpha*(1.0_dp + hpart)/2.0_dp)
614    
615 gezelter 1174 vpair = vpair + pot_temp
616 chuckv 1175
617 gezelter 1174 if (do_pot) then
618     #ifdef IS_MPI
619     pot_row(VDW_POT,atom1) = pot_row(VDW_POT,atom1) + 0.5_dp*pot_temp*sw
620     pot_col(VDW_POT,atom2) = pot_col(VDW_POT,atom2) + 0.5_dp*pot_temp*sw
621     #else
622     pot = pot + pot_temp*sw
623     #endif
624     endif
625    
626     ! derivative wrt r
627 chuckv 1228 dmawdr = 2.0_dp*D0*beta0*(expfnc0 - expfnc02) * &
628     (1.0_dp - alpha*(1.0_dp + hpart)/2.0_dp)
629 gezelter 1174
630     ! derivative wrt cos(theta)
631 chuckv 1228 !dmawdct = 2.0_dp * gamma * D0 * expfncH2 * alpha * sp * sp
632     dmawdct = -sqrt3_4*alpha*D0*vMorse
633    
634     ! derivative wrt sin(theta)
635     ! dmawdst = 2.0_dp * alpha * D0 * vMorse* c2t * st
636     dmawdst = 3.0_dp/4.0_dp*alpha*D0*vmorse*c2p*st
637    
638     ! derivative wrt cos(2*phi)
639     dmawdsp = 3.0_dp/4.0_dp * alpha * D0 * vmorse*st*st
640 gezelter 1174
641 chuckv 1175
642 gezelter 1174 ! chain rule to put these back on x, y, z, ux, uy, uz
643     dvdx = dmawdr * drdx + dmawdct * dctdx + dmawdsp * dspdx
644     dvdy = dmawdr * drdy + dmawdct * dctdy + dmawdsp * dspdy
645     dvdz = dmawdr * drdz + dmawdct * dctdz + dmawdsp * dspdz
646 chuckv 1175
647 gezelter 1174
648 chuckv 1175
649    
650 gezelter 1174 dvdux = dmawdr * drdux + dmawdct * dctdux + dmawdsp * dspdux
651     dvduy = dmawdr * drduy + dmawdct * dctduy + dmawdsp * dspduy
652     dvduz = dmawdr * drduz + dmawdct * dctduz + dmawdsp * dspduz
653    
654 chuckv 1175
655    
656    
657 gezelter 1174 tx = (dvduy - dvduz) * sw
658     ty = (dvduz - dvdux) * sw
659     tz = (dvdux - dvduy) * sw
660    
661     ! go back to lab frame using transpose of rotation matrix:
662    
663     #ifdef IS_MPI
664     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
665     t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*tx + &
666     a_Row(4,atom1)*ty + a_Row(7,atom1)*tz
667     t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*tx + &
668     a_Row(5,atom1)*ty + a_Row(8,atom1)*tz
669     t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*tx + &
670     a_Row(6,atom1)*ty + a_Row(9,atom1)*tz
671     else
672     t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*tx + &
673     a_Col(4,atom2)*ty + a_Col(7,atom2)*tz
674     t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*tx + &
675     a_Col(5,atom2)*ty + a_Col(8,atom2)*tz
676     t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*tx + &
677     a_Col(6,atom2)*ty + a_Col(9,atom2)*tz
678     endif
679     #else
680     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
681     t(1,atom1) = t(1,atom1) + a(1,atom1)*tx + a(4,atom1)*ty + &
682     a(7,atom1)*tz
683     t(2,atom1) = t(2,atom1) + a(2,atom1)*tx + a(5,atom1)*ty + &
684     a(8,atom1)*tz
685     t(3,atom1) = t(3,atom1) + a(3,atom1)*tx + a(6,atom1)*ty + &
686     a(9,atom1)*tz
687     else
688     t(1,atom2) = t(1,atom2) + a(1,atom2)*tx + a(4,atom2)*ty + &
689     a(7,atom2)*tz
690     t(2,atom2) = t(2,atom2) + a(2,atom2)*tx + a(5,atom2)*ty + &
691     a(8,atom2)*tz
692     t(3,atom2) = t(3,atom2) + a(3,atom2)*tx + a(6,atom2)*ty + &
693     a(9,atom2)*tz
694     endif
695     #endif
696     ! Now, on to the forces:
697    
698     fx = -dvdx * sw
699     fy = -dvdy * sw
700     fz = -dvdz * sw
701    
702     ! rotate the terms back into the lab frame:
703    
704     #ifdef IS_MPI
705     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
706     fxl = a_Row(1,atom1)*fx + a_Row(4,atom1)*fy + a_Row(7,atom1)*fz
707     fyl = a_Row(2,atom1)*fx + a_Row(5,atom1)*fy + a_Row(8,atom1)*fz
708     fzl = a_Row(3,atom1)*fx + a_Row(6,atom1)*fy + a_Row(9,atom1)*fz
709     else
710     fxl = a_Col(1,atom2)*fx + a_Col(4,atom2)*fy + a_Col(7,atom2)*fz
711     fyl = a_Col(2,atom2)*fx + a_Col(5,atom2)*fy + a_Col(8,atom2)*fz
712     fzl = a_Col(3,atom2)*fx + a_Col(6,atom2)*fy + a_Col(9,atom2)*fz
713     endif
714     f_Row(1,atom1) = f_Row(1,atom1) + fxl
715     f_Row(2,atom1) = f_Row(2,atom1) + fyl
716     f_Row(3,atom1) = f_Row(3,atom1) + fzl
717    
718     f_Col(1,atom2) = f_Col(1,atom2) - fxl
719     f_Col(2,atom2) = f_Col(2,atom2) - fyl
720     f_Col(3,atom2) = f_Col(3,atom2) - fzl
721     #else
722     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
723     fxl = a(1,atom1)*fx + a(4,atom1)*fy + a(7,atom1)*fz
724     fyl = a(2,atom1)*fx + a(5,atom1)*fy + a(8,atom1)*fz
725     fzl = a(3,atom1)*fx + a(6,atom1)*fy + a(9,atom1)*fz
726     else
727     fxl = a(1,atom2)*fx + a(4,atom2)*fy + a(7,atom2)*fz
728     fyl = a(2,atom2)*fx + a(5,atom2)*fy + a(8,atom2)*fz
729     fzl = a(3,atom2)*fx + a(6,atom2)*fy + a(9,atom2)*fz
730     endif
731     f(1,atom1) = f(1,atom1) + fxl
732     f(2,atom1) = f(2,atom1) + fyl
733     f(3,atom1) = f(3,atom1) + fzl
734    
735     f(1,atom2) = f(1,atom2) - fxl
736     f(2,atom2) = f(2,atom2) - fyl
737     f(3,atom2) = f(3,atom2) - fzl
738     #endif
739    
740     #ifdef IS_MPI
741     id1 = AtomRowToGlobal(atom1)
742     id2 = AtomColToGlobal(atom2)
743     #else
744     id1 = atom1
745     id2 = atom2
746     #endif
747    
748     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
749    
750     fpair(1) = fpair(1) + fxl
751     fpair(2) = fpair(2) + fyl
752     fpair(3) = fpair(3) + fzl
753    
754     endif
755    
756     return
757 chuckv 1173 end subroutine calc_mnm_maw
758    
759    
760 chuckv 1168 subroutine setMnMDefaultCutoff(thisRcut, shiftedPot, shiftedFrc)
761     real(kind=dp), intent(in) :: thisRcut
762     logical, intent(in) :: shiftedPot
763     logical, intent(in) :: shiftedFrc
764     integer i, nInteractions
765     defaultCutoff = thisRcut
766     defaultShiftPot = shiftedPot
767     defaultShiftFrc = shiftedFrc
768    
769 xsun 1178 if (associated(MnM_Map)) then
770     if(MnM_Map%interactionCount /= 0) then
771     nInteractions = MnM_Map%interactionCount
772 chuckv 1168
773 xsun 1178 do i = 1, nInteractions
774     MnM_Map%interactions(i)%shiftedPot = shiftedPot
775     MnM_Map%interactions(i)%shiftedFrc = shiftedFrc
776     MnM_Map%interactions(i)%rCut = thisRcut
777     MnM_Map%interactions(i)%rCutWasSet = .true.
778     enddo
779     end if
780     end if
781 chuckv 1168
782     end subroutine setMnMDefaultCutoff
783    
784     subroutine copyAllData(v1, v2)
785     type(MnMinteractionMap), pointer :: v1
786     type(MnMinteractionMap), pointer :: v2
787     integer :: i, j
788    
789     do i = 1, v1%interactionCount
790     v2%interactions(i) = v1%interactions(i)
791     enddo
792    
793     v2%interactionCount = v1%interactionCount
794     return
795     end subroutine copyAllData
796    
797     subroutine addInteraction(myInteraction)
798     type(MNMtype) :: myInteraction
799     type(MnMinteraction) :: nt
800     integer :: id
801 chuckv 1175
802    
803 chuckv 1168
804 chuckv 1175
805    
806 chuckv 1168 nt%interaction_type = myInteraction%MNMInteractionType
807 chuckv 1175 nt%metal_atid = &
808     getFirstMatchingElement(atypes, "c_ident", myInteraction%metal_atid)
809     nt%nonmetal_atid = &
810     getFirstMatchingElement(atypes, "c_ident", myInteraction%nonmetal_atid)
811 chuckv 1168
812 chuckv 1175
813 chuckv 1168 select case (nt%interaction_type)
814     case (MNM_LENNARDJONES)
815     nt%sigma = myInteraction%sigma
816     nt%epsilon = myInteraction%epsilon
817     case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE)
818     nt%R0 = myInteraction%R0
819     nt%D0 = myInteraction%D0
820     nt%beta0 = myInteraction%beta0
821     case(MNM_MAW)
822     nt%R0 = myInteraction%R0
823     nt%D0 = myInteraction%D0
824     nt%beta0 = myInteraction%beta0
825     nt%betaH = myInteraction%betaH
826     nt%alpha = myInteraction%alpha
827     nt%gamma = myInteraction%gamma
828     case default
829 chuckv 1173 call handleError("MNM", "Unknown Interaction type")
830 chuckv 1168 end select
831    
832     if (.not. associated(MnM_Map)) then
833     call ensureCapacityHelper(MnM_Map, 1)
834     else
835     call ensureCapacityHelper(MnM_Map, MnM_Map%interactionCount + 1)
836     end if
837    
838     MnM_Map%interactionCount = MnM_Map%interactionCount + 1
839     id = MnM_Map%interactionCount
840     MnM_Map%interactions(id) = nt
841     end subroutine addInteraction
842    
843     subroutine ensureCapacityHelper(this, minCapacity)
844     type(MnMinteractionMap), pointer :: this, that
845     integer, intent(in) :: minCapacity
846     integer :: oldCapacity
847     integer :: newCapacity
848     logical :: resizeFlag
849    
850     resizeFlag = .false.
851    
852     ! first time: allocate a new vector with default size
853    
854     if (.not. associated(this)) then
855     this => MnMinitialize(minCapacity, 0)
856     endif
857    
858     oldCapacity = size(this%interactions)
859    
860     if (minCapacity > oldCapacity) then
861     if (this%capacityIncrement .gt. 0) then
862     newCapacity = oldCapacity + this%capacityIncrement
863     else
864     newCapacity = oldCapacity * 2
865     endif
866     if (newCapacity .lt. minCapacity) then
867     newCapacity = minCapacity
868     endif
869     resizeFlag = .true.
870     else
871     newCapacity = oldCapacity
872     endif
873    
874     if (resizeFlag) then
875     that => MnMinitialize(newCapacity, this%capacityIncrement)
876     call copyAllData(this, that)
877     this => MnMdestroy(this)
878     this => that
879     endif
880     end subroutine ensureCapacityHelper
881    
882     function MnMinitialize(cap, capinc) result(this)
883     integer, intent(in) :: cap, capinc
884     integer :: error
885     type(MnMinteractionMap), pointer :: this
886    
887     nullify(this)
888    
889     if (cap < 0) then
890     write(*,*) 'Bogus Capacity:', cap
891     return
892     endif
893     allocate(this,stat=error)
894     if ( error /= 0 ) then
895     write(*,*) 'Could not allocate MnMinteractionMap!'
896     return
897     end if
898    
899     this%initialCapacity = cap
900     this%capacityIncrement = capinc
901    
902     allocate(this%interactions(this%initialCapacity), stat=error)
903     if(error /= 0) write(*,*) 'Could not allocate MnMinteraction!'
904    
905     end function MnMinitialize
906    
907 chuckv 1175 subroutine createInteractionLookup(this)
908     type (MnMInteractionMap),pointer :: this
909 chuckv 1168 integer :: biggestAtid, i, metal_atid, nonmetal_atid, error
910    
911     biggestAtid=-1
912     do i = 1, this%interactionCount
913     metal_atid = this%interactions(i)%metal_atid
914     nonmetal_atid = this%interactions(i)%nonmetal_atid
915    
916     if (metal_atid .gt. biggestAtid) biggestAtid = metal_atid
917     if (nonmetal_atid .gt. biggestAtid) biggestAtid = nonmetal_atid
918     enddo
919 chuckv 1175
920 chuckv 1168 allocate(MnMinteractionLookup(biggestAtid,biggestAtid), stat=error)
921     if (error /= 0) write(*,*) 'Could not allocate MnMinteractionLookup'
922    
923     do i = 1, this%interactionCount
924     metal_atid = this%interactions(i)%metal_atid
925     nonmetal_atid = this%interactions(i)%nonmetal_atid
926    
927     MnMinteractionLookup(metal_atid, nonmetal_atid) = i
928     MnMinteractionLookup(nonmetal_atid, metal_atid) = i
929     enddo
930     end subroutine createInteractionLookup
931    
932    
933     function MnMdestroy(this) result(null_this)
934     logical :: done
935     type(MnMinteractionMap), pointer :: this
936     type(MnMinteractionMap), pointer :: null_this
937    
938     if (.not. associated(this)) then
939     null_this => null()
940     return
941     end if
942    
943     !! Walk down the list and deallocate each of the map's components
944     if(associated(this%interactions)) then
945     deallocate(this%interactions)
946     this%interactions=>null()
947     endif
948     deallocate(this)
949     this => null()
950     null_this => null()
951     end function MnMdestroy
952    
953    
954     subroutine deleteInteractions()
955     MnM_Map => MnMdestroy(MnM_Map)
956     return
957     end subroutine deleteInteractions
958    
959 chuckv 1173
960     subroutine getLJfunc(r, myPot, myDeriv)
961    
962     real(kind=dp), intent(in) :: r
963     real(kind=dp), intent(inout) :: myPot, myDeriv
964     real(kind=dp) :: ri, ri2, ri6, ri7, ri12, ri13
965     real(kind=dp) :: a, b, c, d, dx
966     integer :: j
967    
968     ri = 1.0_DP / r
969     ri2 = ri*ri
970     ri6 = ri2*ri2*ri2
971     ri7 = ri6*ri
972     ri12 = ri6*ri6
973     ri13 = ri12*ri
974    
975     myPot = 4.0_DP * (ri12 - ri6)
976     myDeriv = 24.0_DP * (ri7 - 2.0_DP * ri13)
977    
978     return
979     end subroutine getLJfunc
980    
981     subroutine getSoftFunc(r, myPot, myDeriv)
982    
983     real(kind=dp), intent(in) :: r
984     real(kind=dp), intent(inout) :: myPot, myDeriv
985     real(kind=dp) :: ri, ri2, ri6, ri7
986     real(kind=dp) :: a, b, c, d, dx
987     integer :: j
988    
989     ri = 1.0_DP / r
990     ri2 = ri*ri
991     ri6 = ri2*ri2*ri2
992     ri7 = ri6*ri
993     myPot = 4.0_DP * (ri6)
994     myDeriv = - 24.0_DP * ri7
995    
996     return
997     end subroutine getSoftFunc
998    
999    
1000    
1001    
1002    
1003    
1004 chuckv 1168 end module MetalNonMetal