ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 1174
Committed: Tue Jul 17 21:55:57 2007 UTC (18 years ago) by gezelter
Original Path: trunk/src/UseTheForce/DarkSide/MetalNonMetal.F90
File size: 29102 byte(s)
Log Message:
adding MetalNonMetal interactions

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 gezelter 1174 !! @version $Id: MetalNonMetal.F90,v 1.3 2007-07-17 21:55:57 gezelter Exp $, $Date: 2007-07-17 21:55:57 $, $Name: not supported by cvs2svn $, $Revision: 1.3 $
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     logical, save :: defaultShiftPot = .false.
71     logical, save :: defaultShiftFrc = .false.
72     logical, save :: haveDefaultCutoff = .false.
73    
74     type :: MnMinteraction
75     integer :: metal_atid
76     integer :: nonmetal_atid
77     integer :: interaction_type
78     real(kind=dp) :: R0
79     real(kind=dp) :: D0
80     real(kind=dp) :: beta0
81     real(kind=dp) :: betaH
82     real(kind=dp) :: alpha
83     real(kind=dp) :: gamma
84     real(kind=dp) :: sigma
85     real(kind=dp) :: epsilon
86     real(kind=dp) :: rCut = 0.0_dp
87     logical :: rCutWasSet = .false.
88     logical :: shiftedPot = .false.
89     logical :: shiftedFrc = .false.
90     end type MnMinteraction
91    
92     type :: MnMinteractionMap
93     PRIVATE
94     integer :: initialCapacity = 10
95     integer :: capacityIncrement = 0
96     integer :: interactionCount = 0
97     type(MnMinteraction), pointer :: interactions(:) => null()
98     end type MnMinteractionMap
99    
100     type (MnMInteractionMap), pointer :: MnM_Map
101    
102     integer, allocatable, dimension(:,:) :: MnMinteractionLookup
103    
104     public :: setMnMDefaultCutoff
105     public :: addInteraction
106     public :: deleteInteractions
107     public :: MNMtype
108 chuckv 1173 public :: do_mnm_pair
109 chuckv 1168
110     contains
111    
112    
113 chuckv 1173 subroutine do_mnm_pair(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
114     Pot, A, F,t, Do_pot)
115 gezelter 1174 integer, intent(in) :: atom1, atom2
116 chuckv 1173 integer :: atid1, atid2, ljt1, ljt2
117 gezelter 1174 real( kind = dp ), intent(in) :: rij, r2, rcut
118     real( kind = dp ) :: pot, sw, vpair
119 chuckv 1173 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
120 gezelter 1174 real (kind=dp), intent(in), dimension(9,nLocal) :: A
121     real (kind=dp), intent(inout), dimension(3,nLocal) :: t
122     real( kind = dp ), intent(in), dimension(3) :: d
123 chuckv 1173 real( kind = dp ), intent(inout), dimension(3) :: fpair
124 gezelter 1174 logical, intent(in) :: do_pot
125 chuckv 1168
126 chuckv 1173 integer :: interaction_id
127     integer :: interaction_type
128    
129     #ifdef IS_MPI
130     atid1 = atid_Row(atom1)
131     atid2 = atid_Col(atom2)
132     #else
133     atid1 = atid(atom1)
134     atid2 = atid(atom2)
135     #endif
136    
137     interaction_id = MnMinteractionLookup(atid1, atid2)
138     interaction_type = MnM_Map%interactions(interaction_id)%interaction_type
139    
140     select case (interaction_type)
141     case (MNM_LENNARDJONES)
142 gezelter 1174 call calc_mnm_lennardjones(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, &
143     Fpair, Pot, F, Do_pot, interaction_id)
144 chuckv 1173 case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE)
145     call calc_mnm_morse(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
146 gezelter 1174 Pot, F, Do_pot, interaction_id, interaction_type)
147 chuckv 1173 case(MNM_MAW)
148     call calc_mnm_maw(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
149 gezelter 1174 Pot,A, F,t, Do_pot, interaction_id)
150 chuckv 1173 end select
151    
152     end subroutine do_mnm_pair
153    
154 gezelter 1174 subroutine calc_mnm_lennardjones(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, &
155     Fpair, Pot, F, Do_pot, interaction_id)
156 chuckv 1173
157 gezelter 1174 integer, intent(in) :: atom1, atom2
158     real( kind = dp ), intent(in) :: rij, r2, rcut
159     real( kind = dp ) :: pot, sw, vpair
160 chuckv 1173 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
161 gezelter 1174 real( kind = dp ), intent(in), dimension(3) :: d
162 chuckv 1173 real( kind = dp ), intent(inout), dimension(3) :: fpair
163 gezelter 1174 logical, intent(in) :: do_pot
164 chuckv 1173 integer, intent(in) :: interaction_id
165    
166     ! local Variables
167     real( kind = dp ) :: drdx, drdy, drdz
168     real( kind = dp ) :: fx, fy, fz
169     real( kind = dp ) :: myPot, myPotC, myDeriv, myDerivC, ros, rcos
170     real( kind = dp ) :: pot_temp, dudr
171     real( kind = dp ) :: sigmai
172     real( kind = dp ) :: epsilon
173     logical :: isSoftCore, shiftedPot, shiftedFrc
174     integer :: id1, id2, localError
175    
176 gezelter 1174 sigmai = 1.0_dp / MnM_Map%interactions(interaction_id)%sigma
177 chuckv 1173 epsilon = MnM_Map%interactions(interaction_id)%epsilon
178     shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot
179     shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc
180    
181     ros = rij * sigmai
182    
183     call getLJfunc(ros, myPot, myDeriv)
184    
185     if (shiftedPot) then
186     rcos = rcut * sigmai
187     call getLJfunc(rcos, myPotC, myDerivC)
188     myDerivC = 0.0_dp
189     elseif (shiftedFrc) then
190     rcos = rcut * sigmai
191     call getLJfunc(rcos, myPotC, myDerivC)
192     myPotC = myPotC + myDerivC * (rij - rcut) * sigmai
193     else
194     myPotC = 0.0_dp
195     myDerivC = 0.0_dp
196 gezelter 1174 endif
197 chuckv 1173
198     pot_temp = epsilon * (myPot - myPotC)
199     vpair = vpair + pot_temp
200     dudr = sw * epsilon * (myDeriv - myDerivC) * sigmai
201    
202     drdx = d(1) / rij
203     drdy = d(2) / rij
204     drdz = d(3) / rij
205    
206     fx = dudr * drdx
207     fy = dudr * drdy
208     fz = dudr * drdz
209 gezelter 1174
210 chuckv 1173 #ifdef IS_MPI
211     if (do_pot) then
212     pot_Row(VDW_POT,atom1) = pot_Row(VDW_POT,atom1) + sw*pot_temp*0.5
213     pot_Col(VDW_POT,atom2) = pot_Col(VDW_POT,atom2) + sw*pot_temp*0.5
214     endif
215 gezelter 1174
216 chuckv 1173 f_Row(1,atom1) = f_Row(1,atom1) + fx
217     f_Row(2,atom1) = f_Row(2,atom1) + fy
218     f_Row(3,atom1) = f_Row(3,atom1) + fz
219 gezelter 1174
220 chuckv 1173 f_Col(1,atom2) = f_Col(1,atom2) - fx
221     f_Col(2,atom2) = f_Col(2,atom2) - fy
222     f_Col(3,atom2) = f_Col(3,atom2) - fz
223 gezelter 1174
224 chuckv 1173 #else
225     if (do_pot) pot = pot + sw*pot_temp
226 gezelter 1174
227 chuckv 1173 f(1,atom1) = f(1,atom1) + fx
228     f(2,atom1) = f(2,atom1) + fy
229     f(3,atom1) = f(3,atom1) + fz
230 gezelter 1174
231 chuckv 1173 f(1,atom2) = f(1,atom2) - fx
232     f(2,atom2) = f(2,atom2) - fy
233     f(3,atom2) = f(3,atom2) - fz
234     #endif
235    
236     #ifdef IS_MPI
237     id1 = AtomRowToGlobal(atom1)
238     id2 = AtomColToGlobal(atom2)
239     #else
240     id1 = atom1
241     id2 = atom2
242     #endif
243 gezelter 1174
244 chuckv 1173 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
245 gezelter 1174
246 chuckv 1173 fpair(1) = fpair(1) + fx
247     fpair(2) = fpair(2) + fy
248     fpair(3) = fpair(3) + fz
249 gezelter 1174
250 chuckv 1173 endif
251 gezelter 1174 return
252 chuckv 1173 end subroutine calc_mnm_lennardjones
253    
254     subroutine calc_mnm_morse(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
255     Pot, f, Do_pot, interaction_id, interaction_type)
256 gezelter 1174 integer, intent(in) :: atom1, atom2
257     real( kind = dp ), intent(in) :: rij, r2, rcut
258     real( kind = dp ) :: pot, sw, vpair
259 chuckv 1173 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
260 gezelter 1174 real( kind = dp ), intent(in), dimension(3) :: d
261 chuckv 1173 real( kind = dp ), intent(inout), dimension(3) :: fpair
262 gezelter 1174 logical, intent(in) :: do_pot
263 chuckv 1173 integer, intent(in) :: interaction_id, interaction_type
264 gezelter 1174 logical :: shiftedPot, shiftedFrc
265 chuckv 1173
266     ! Local Variables
267     real(kind=dp) :: Beta0
268     real(kind=dp) :: R0
269     real(kind=dp) :: D0
270     real(kind=dp) :: expt
271     real(kind=dp) :: expt2
272     real(kind=dp) :: expfnc
273     real(kind=dp) :: expfnc2
274     real(kind=dp) :: D_expt
275     real(kind=dp) :: D_expt2
276     real(kind=dp) :: rcos
277     real(kind=dp) :: exptC
278     real(kind=dp) :: expt2C
279     real(kind=dp) :: expfncC
280     real(kind=dp) :: expfnc2C
281     real(kind=dp) :: D_expt2C
282     real(kind=dp) :: D_exptC
283    
284     real(kind=dp) :: myPot
285     real(kind=dp) :: myPotC
286     real(kind=dp) :: myDeriv
287     real(kind=dp) :: myDerivC
288     real(kind=dp) :: pot_temp
289     real(kind=dp) :: fx,fy,fz
290     real(kind=dp) :: drdx,drdy,drdz
291     real(kind=dp) :: dudr
292     integer :: id1,id2
293 gezelter 1174
294 chuckv 1173
295 gezelter 1174 D0 = MnM_Map%interactions(interaction_id)%D0
296     R0 = MnM_Map%interactions(interaction_id)%r0
297 chuckv 1173 Beta0 = MnM_Map%interactions(interaction_id)%Beta0
298 gezelter 1174 shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot
299     shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc
300    
301     ! V(r) = D_e exp(-a(r-re)(exp(-a(r-re))-2)
302    
303     expt = -beta0*(rij - R0)
304 chuckv 1173 expfnc = exp(expt)
305 gezelter 1174 expfnc2 = expfnc*expfnc
306 chuckv 1173
307 gezelter 1174 if (shiftedPot .or. shiftedFrc) then
308     exptC = -beta0*(rcut - R0)
309     expfncC = exp(exptC)
310     expfnc2C = expfncC*expfncC
311     endif
312    
313 chuckv 1173 select case (interaction_type)
314 gezelter 1174 case (MNM_SHIFTEDMORSE)
315 chuckv 1173
316 gezelter 1174 myPot = D0 * (expfnc2 - 2.0_dp * expfnc)
317     myDeriv = 2.0*D0*beta0*(expfnc - expfnc2)
318 chuckv 1173
319 gezelter 1174 if (shiftedPot) then
320     myPotC = D0 * (expfnc2C - 2.0_dp * expfncC)
321     myDerivC = 0.0_dp
322     elseif (shiftedFrc) then
323     myPotC = D0 * (expfnc2C - 2.0_dp * expfncC)
324     myDerivC = 2.0*D0*beta0*(expfnc2C - expfnc2C)
325     myPotC = myPotC + myDerivC * (rij - rcut)
326     else
327     myPotC = 0.0_dp
328     myDerivC = 0.0_dp
329     endif
330 chuckv 1173
331     case (MNM_REPULSIVEMORSE)
332    
333 gezelter 1174 myPot = D0 * expfnc2
334     myDeriv = -2.0_dp * D0 * beta0 * expfnc2
335 chuckv 1173
336 gezelter 1174 if (shiftedPot) then
337     myPotC = D0 * expfnc2C
338     myDerivC = 0.0_dp
339     elseif (shiftedFrc) then
340     myPotC = D0 * expfnc2C
341     myDerivC = -2.0_dp * D0* beta0 * expfnc2C
342     myPotC = myPotC + myDerivC * (rij - rcut)
343     else
344     myPotC = 0.0_dp
345     myDerivC = 0.0_dp
346     endif
347 chuckv 1173 end select
348    
349     pot_temp = (myPot - myPotC)
350     vpair = vpair + pot_temp
351     dudr = sw * (myDeriv - myDerivC)
352    
353     drdx = d(1) / rij
354     drdy = d(2) / rij
355     drdz = d(3) / rij
356    
357     fx = dudr * drdx
358     fy = dudr * drdy
359     fz = dudr * drdz
360    
361     #ifdef IS_MPI
362     if (do_pot) then
363     pot_Row(VDW_POT,atom1) = pot_Row(VDW_POT,atom1) + sw*pot_temp*0.5
364     pot_Col(VDW_POT,atom2) = pot_Col(VDW_POT,atom2) + sw*pot_temp*0.5
365     endif
366    
367     f_Row(1,atom1) = f_Row(1,atom1) + fx
368     f_Row(2,atom1) = f_Row(2,atom1) + fy
369     f_Row(3,atom1) = f_Row(3,atom1) + fz
370    
371     f_Col(1,atom2) = f_Col(1,atom2) - fx
372     f_Col(2,atom2) = f_Col(2,atom2) - fy
373     f_Col(3,atom2) = f_Col(3,atom2) - fz
374    
375     #else
376     if (do_pot) pot = pot + sw*pot_temp
377    
378     f(1,atom1) = f(1,atom1) + fx
379     f(2,atom1) = f(2,atom1) + fy
380     f(3,atom1) = f(3,atom1) + fz
381    
382     f(1,atom2) = f(1,atom2) - fx
383     f(2,atom2) = f(2,atom2) - fy
384     f(3,atom2) = f(3,atom2) - fz
385     #endif
386    
387     #ifdef IS_MPI
388     id1 = AtomRowToGlobal(atom1)
389     id2 = AtomColToGlobal(atom2)
390     #else
391     id1 = atom1
392     id2 = atom2
393     #endif
394    
395     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
396    
397     fpair(1) = fpair(1) + fx
398     fpair(2) = fpair(2) + fy
399     fpair(3) = fpair(3) + fz
400    
401     endif
402    
403     return
404     end subroutine calc_mnm_morse
405    
406     subroutine calc_mnm_maw(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
407 gezelter 1174 Pot, A, F,t, Do_pot, interaction_id)
408     integer, intent(in) :: atom1, atom2
409     real( kind = dp ), intent(in) :: rij, r2, rcut
410     real( kind = dp ) :: pot, sw, vpair
411 chuckv 1173 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
412 gezelter 1174 real (kind=dp),intent(in), dimension(9,nLocal) :: A
413 chuckv 1173 real (kind=dp),intent(inout), dimension(3,nLocal) :: t
414    
415 gezelter 1174 real( kind = dp ), intent(in), dimension(3) :: d
416 chuckv 1173 real( kind = dp ), intent(inout), dimension(3) :: fpair
417 gezelter 1174 logical, intent(in) :: do_pot
418 chuckv 1173
419     integer, intent(in) :: interaction_id
420    
421 gezelter 1174 real(kind = dp) :: D0, R0, beta0, betaH, alpha, gamma, pot_temp
422     real(kind = dp) :: expt0, expfnc0, expfnc02
423     real(kind = dp) :: exptH, expfncH, expfncH2
424     real(kind = dp) :: x, y, z, x2, y2, z2, r3, proj, proj3, st2
425     real(kind = dp) :: drdx, drdy, drdz, drdux, drduy, drduz
426     real(kind = dp) :: ct, dctdx, dctdy, dctdz, dctdux, dctduy, dctduz
427     real(kind = dp) :: sp, dspdx, dspdy, dspdz, dspdux, dspduy, dspduz
428     real(kind = dp) :: dvdx, dvdy, dvdz, dvdux, dvduy, dvduz
429     real(kind = dp) :: maw, dmawdr, dmawdct, dmawdsp
430     real(kind = dp) :: fx, fy, fz, tx, ty, tz, fxl, fyl, fzl
431     integer :: atid1, atid2, id1, id2
432     logical :: shiftedPot, shiftedFrc
433    
434     #ifdef IS_MPI
435     atid1 = atid_Row(atom1)
436     atid2 = atid_Col(atom2)
437    
438     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
439     ! rotate the inter-particle separation into the two different
440     ! body-fixed coordinate systems:
441    
442     x = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3)
443     y = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3)
444     z = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3)
445     else
446     ! negative sign because this is the vector from j to i:
447    
448     x = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3))
449     y = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3))
450     z = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3))
451     endif
452    
453     #else
454     atid1 = atid(atom1)
455     atid2 = atid(atom2)
456    
457     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
458     ! rotate the inter-particle separation into the two different
459     ! body-fixed coordinate systems:
460    
461     x = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3)
462     y = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3)
463     z = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3)
464     else
465     ! negative sign because this is the vector from j to i:
466    
467     x = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3))
468     y = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3))
469     z = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3))
470     endif
471    
472     #endif
473    
474     D0 = MnM_Map%interactions(interaction_id)%D0
475     R0 = MnM_Map%interactions(interaction_id)%r0
476     beta0 = MnM_Map%interactions(interaction_id)%beta0
477     betaH = MnM_Map%interactions(interaction_id)%betaH
478     alpha = MnM_Map%interactions(interaction_id)%alpha
479     gamma = MnM_Map%interactions(interaction_id)%gamma
480    
481     shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot
482     shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc
483    
484     expt0 = -beta0*(rij - R0)
485     expfnc0 = exp(expt0)
486     expfnc02 = expfnc0*expfnc0
487    
488     exptH = -betaH*(rij - R0)
489     expfncH = exp(exptH)
490     expfncH2 = expfncH*expfncH
491    
492     !!$ if (shiftedPot .or. shiftedFrc) then
493     !!$ exptC0 = -beta0*(rcut - R0)
494     !!$ expfncC0 = exp(exptC0)
495     !!$ expfncC02 = expfncC0*expfncC0
496     !!$ exptCH = -betaH*(rcut - R0)
497     !!$ expfncCH = exp(exptCH)
498     !!$ expfncCH2 = expfncCH*expfncCH
499     !!$ endif
500    
501     drdx = -d(1) / rij
502     drdy = -d(2) / rij
503     drdz = -d(3) / rij
504     drdux = 0.0_dp
505     drduy = 0.0_dp
506     drduz = 0.0_dp
507    
508     x2 = x*x
509     y2 = y*y
510     z2 = z*z
511    
512     ct = z / rij
513    
514     if (ct .gt. 1.0_dp) ct = 1.0_dp
515     if (ct .lt. -1.0_dp) ct = -1.0_dp
516    
517     dctdx = - z * x / r3
518     dctdy = - z * y / r3
519     dctdz = 1.0_dp / rij - z2 / r3
520     dctdux = y / rij ! - (z * x2) / r3
521     dctduy = -x / rij !- (z * y2) / r3
522     dctduz = 0.0_dp ! z / rij - (z2 * z) / r3
523    
524     ! this is an attempt to try to truncate the singularity when
525     ! sin(theta) is near 0.0:
526    
527     st2 = 1.0_dp - ct*ct
528     if (abs(st2) .lt. 1.0e-12_dp) then
529     proj = sqrt(rij * 1.0e-12_dp)
530     !!$ dcpdx = 1.0_dp / proj
531     !!$ dcpdy = 0.0_dp
532     !!$ dcpdux = x / proj
533     !!$ dcpduy = 0.0_dp
534     dspdx = 0.0_dp
535     dspdy = 1.0_dp / proj
536     dspdux = 0.0_dp
537     dspduy = y / proj
538     else
539     proj = sqrt(x2 + y2)
540     proj3 = proj*proj*proj
541     !!$ dcpdx = 1.0_dp / proj - x2 / proj3
542     !!$ dcpdy = - x * y / proj3
543     !!$ dcpdux = x / proj - (x2 * x) / proj3
544     !!$ dcpduy = - (x * y2) / proj3
545     dspdx = - x * y / proj3
546     dspdy = 1.0_dp / proj - y2 / proj3
547     dspdux = - (y * x2) / proj3
548     dspduy = y / proj - (y2 * y) / proj3
549     endif
550    
551     !!$ cp = x / proj
552     !!$ dcpdz = 0.0_dp
553     !!$ dcpduz = 0.0_dp
554    
555     sp = y / proj
556     dspdz = 0.0_dp
557     dspduz = 0.0_dp
558    
559    
560     pot_temp = D0 * (expfnc02 - 2.0_dp * expfnc0) + &
561     2.0_dp * gamma * D0 * expfncH2 * (1.0_dp + alpha*ct)*sp*sp
562     vpair = vpair + pot_temp
563     if (do_pot) then
564     #ifdef IS_MPI
565     pot_row(VDW_POT,atom1) = pot_row(VDW_POT,atom1) + 0.5_dp*pot_temp*sw
566     pot_col(VDW_POT,atom2) = pot_col(VDW_POT,atom2) + 0.5_dp*pot_temp*sw
567     #else
568     pot = pot + pot_temp*sw
569     #endif
570     endif
571    
572     ! derivative wrt r
573     dmawdr = 2.0*D0*beta0*(expfnc0 - expfnc02) - &
574     4.0 * gamma * D0 * betaH * expfncH2 * (1.0_dp + alpha*ct)*sp*sp
575    
576     ! derivative wrt cos(theta)
577     dmawdct = 2.0 * gamma * D0 * expfncH2 * alpha * sp * sp
578    
579     ! derivative wrt sin(phi)
580     dmawdsp = 4.0 * gamma * D0 * expfncH2 * (1.0_dp + alpha*ct)*sp
581    
582     ! chain rule to put these back on x, y, z, ux, uy, uz
583     dvdx = dmawdr * drdx + dmawdct * dctdx + dmawdsp * dspdx
584     dvdy = dmawdr * drdy + dmawdct * dctdy + dmawdsp * dspdy
585     dvdz = dmawdr * drdz + dmawdct * dctdz + dmawdsp * dspdz
586    
587     dvdux = dmawdr * drdux + dmawdct * dctdux + dmawdsp * dspdux
588     dvduy = dmawdr * drduy + dmawdct * dctduy + dmawdsp * dspduy
589     dvduz = dmawdr * drduz + dmawdct * dctduz + dmawdsp * dspduz
590    
591     tx = (dvduy - dvduz) * sw
592     ty = (dvduz - dvdux) * sw
593     tz = (dvdux - dvduy) * sw
594    
595     ! go back to lab frame using transpose of rotation matrix:
596    
597     #ifdef IS_MPI
598     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
599     t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*tx + &
600     a_Row(4,atom1)*ty + a_Row(7,atom1)*tz
601     t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*tx + &
602     a_Row(5,atom1)*ty + a_Row(8,atom1)*tz
603     t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*tx + &
604     a_Row(6,atom1)*ty + a_Row(9,atom1)*tz
605     else
606     t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*tx + &
607     a_Col(4,atom2)*ty + a_Col(7,atom2)*tz
608     t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*tx + &
609     a_Col(5,atom2)*ty + a_Col(8,atom2)*tz
610     t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*tx + &
611     a_Col(6,atom2)*ty + a_Col(9,atom2)*tz
612     endif
613     #else
614     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
615     t(1,atom1) = t(1,atom1) + a(1,atom1)*tx + a(4,atom1)*ty + &
616     a(7,atom1)*tz
617     t(2,atom1) = t(2,atom1) + a(2,atom1)*tx + a(5,atom1)*ty + &
618     a(8,atom1)*tz
619     t(3,atom1) = t(3,atom1) + a(3,atom1)*tx + a(6,atom1)*ty + &
620     a(9,atom1)*tz
621     else
622     t(1,atom2) = t(1,atom2) + a(1,atom2)*tx + a(4,atom2)*ty + &
623     a(7,atom2)*tz
624     t(2,atom2) = t(2,atom2) + a(2,atom2)*tx + a(5,atom2)*ty + &
625     a(8,atom2)*tz
626     t(3,atom2) = t(3,atom2) + a(3,atom2)*tx + a(6,atom2)*ty + &
627     a(9,atom2)*tz
628     endif
629     #endif
630     ! Now, on to the forces:
631    
632     fx = -dvdx * sw
633     fy = -dvdy * sw
634     fz = -dvdz * sw
635    
636     ! rotate the terms back into the lab frame:
637    
638     #ifdef IS_MPI
639     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
640     fxl = a_Row(1,atom1)*fx + a_Row(4,atom1)*fy + a_Row(7,atom1)*fz
641     fyl = a_Row(2,atom1)*fx + a_Row(5,atom1)*fy + a_Row(8,atom1)*fz
642     fzl = a_Row(3,atom1)*fx + a_Row(6,atom1)*fy + a_Row(9,atom1)*fz
643     else
644     fxl = a_Col(1,atom2)*fx + a_Col(4,atom2)*fy + a_Col(7,atom2)*fz
645     fyl = a_Col(2,atom2)*fx + a_Col(5,atom2)*fy + a_Col(8,atom2)*fz
646     fzl = a_Col(3,atom2)*fx + a_Col(6,atom2)*fy + a_Col(9,atom2)*fz
647     endif
648     f_Row(1,atom1) = f_Row(1,atom1) + fxl
649     f_Row(2,atom1) = f_Row(2,atom1) + fyl
650     f_Row(3,atom1) = f_Row(3,atom1) + fzl
651    
652     f_Col(1,atom2) = f_Col(1,atom2) - fxl
653     f_Col(2,atom2) = f_Col(2,atom2) - fyl
654     f_Col(3,atom2) = f_Col(3,atom2) - fzl
655     #else
656     if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
657     fxl = a(1,atom1)*fx + a(4,atom1)*fy + a(7,atom1)*fz
658     fyl = a(2,atom1)*fx + a(5,atom1)*fy + a(8,atom1)*fz
659     fzl = a(3,atom1)*fx + a(6,atom1)*fy + a(9,atom1)*fz
660     else
661     fxl = a(1,atom2)*fx + a(4,atom2)*fy + a(7,atom2)*fz
662     fyl = a(2,atom2)*fx + a(5,atom2)*fy + a(8,atom2)*fz
663     fzl = a(3,atom2)*fx + a(6,atom2)*fy + a(9,atom2)*fz
664     endif
665     f(1,atom1) = f(1,atom1) + fxl
666     f(2,atom1) = f(2,atom1) + fyl
667     f(3,atom1) = f(3,atom1) + fzl
668    
669     f(1,atom2) = f(1,atom2) - fxl
670     f(2,atom2) = f(2,atom2) - fyl
671     f(3,atom2) = f(3,atom2) - fzl
672     #endif
673    
674     #ifdef IS_MPI
675     id1 = AtomRowToGlobal(atom1)
676     id2 = AtomColToGlobal(atom2)
677     #else
678     id1 = atom1
679     id2 = atom2
680     #endif
681    
682     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
683    
684     fpair(1) = fpair(1) + fxl
685     fpair(2) = fpair(2) + fyl
686     fpair(3) = fpair(3) + fzl
687    
688     endif
689    
690     return
691 chuckv 1173 end subroutine calc_mnm_maw
692    
693    
694 chuckv 1168 subroutine setMnMDefaultCutoff(thisRcut, shiftedPot, shiftedFrc)
695     real(kind=dp), intent(in) :: thisRcut
696     logical, intent(in) :: shiftedPot
697     logical, intent(in) :: shiftedFrc
698     integer i, nInteractions
699     defaultCutoff = thisRcut
700     defaultShiftPot = shiftedPot
701     defaultShiftFrc = shiftedFrc
702    
703     if(MnM_Map%interactionCount /= 0) then
704     nInteractions = MnM_Map%interactionCount
705    
706     do i = 1, nInteractions
707     MnM_Map%interactions(i)%shiftedPot = shiftedPot
708     MnM_Map%interactions(i)%shiftedFrc = shiftedFrc
709     MnM_Map%interactions(i)%rCut = thisRcut
710     MnM_Map%interactions(i)%rCutWasSet = .true.
711     enddo
712     end if
713    
714     end subroutine setMnMDefaultCutoff
715    
716     subroutine copyAllData(v1, v2)
717     type(MnMinteractionMap), pointer :: v1
718     type(MnMinteractionMap), pointer :: v2
719     integer :: i, j
720    
721     do i = 1, v1%interactionCount
722     v2%interactions(i) = v1%interactions(i)
723     enddo
724    
725     v2%interactionCount = v1%interactionCount
726     return
727     end subroutine copyAllData
728    
729     subroutine addInteraction(myInteraction)
730     type(MNMtype) :: myInteraction
731     type(MnMinteraction) :: nt
732     integer :: id
733    
734     nt%interaction_type = myInteraction%MNMInteractionType
735     nt%metal_atid = myInteraction%metal_atid
736     nt%nonmetal_atid = myInteraction%nonmetal_atid
737    
738     select case (nt%interaction_type)
739     case (MNM_LENNARDJONES)
740     nt%sigma = myInteraction%sigma
741     nt%epsilon = myInteraction%epsilon
742     case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE)
743     nt%R0 = myInteraction%R0
744     nt%D0 = myInteraction%D0
745     nt%beta0 = myInteraction%beta0
746     case(MNM_MAW)
747     nt%R0 = myInteraction%R0
748     nt%D0 = myInteraction%D0
749     nt%beta0 = myInteraction%beta0
750     nt%betaH = myInteraction%betaH
751     nt%alpha = myInteraction%alpha
752     nt%gamma = myInteraction%gamma
753     case default
754 chuckv 1173 call handleError("MNM", "Unknown Interaction type")
755 chuckv 1168 end select
756    
757     if (.not. associated(MnM_Map)) then
758     call ensureCapacityHelper(MnM_Map, 1)
759     else
760     call ensureCapacityHelper(MnM_Map, MnM_Map%interactionCount + 1)
761     end if
762    
763     MnM_Map%interactionCount = MnM_Map%interactionCount + 1
764     id = MnM_Map%interactionCount
765     MnM_Map%interactions(id) = nt
766     end subroutine addInteraction
767    
768     subroutine ensureCapacityHelper(this, minCapacity)
769     type(MnMinteractionMap), pointer :: this, that
770     integer, intent(in) :: minCapacity
771     integer :: oldCapacity
772     integer :: newCapacity
773     logical :: resizeFlag
774    
775     resizeFlag = .false.
776    
777     ! first time: allocate a new vector with default size
778    
779     if (.not. associated(this)) then
780     this => MnMinitialize(minCapacity, 0)
781     endif
782    
783     oldCapacity = size(this%interactions)
784    
785     if (minCapacity > oldCapacity) then
786     if (this%capacityIncrement .gt. 0) then
787     newCapacity = oldCapacity + this%capacityIncrement
788     else
789     newCapacity = oldCapacity * 2
790     endif
791     if (newCapacity .lt. minCapacity) then
792     newCapacity = minCapacity
793     endif
794     resizeFlag = .true.
795     else
796     newCapacity = oldCapacity
797     endif
798    
799     if (resizeFlag) then
800     that => MnMinitialize(newCapacity, this%capacityIncrement)
801     call copyAllData(this, that)
802     this => MnMdestroy(this)
803     this => that
804     endif
805     end subroutine ensureCapacityHelper
806    
807     function MnMinitialize(cap, capinc) result(this)
808     integer, intent(in) :: cap, capinc
809     integer :: error
810     type(MnMinteractionMap), pointer :: this
811    
812     nullify(this)
813    
814     if (cap < 0) then
815     write(*,*) 'Bogus Capacity:', cap
816     return
817     endif
818     allocate(this,stat=error)
819     if ( error /= 0 ) then
820     write(*,*) 'Could not allocate MnMinteractionMap!'
821     return
822     end if
823    
824     this%initialCapacity = cap
825     this%capacityIncrement = capinc
826    
827     allocate(this%interactions(this%initialCapacity), stat=error)
828     if(error /= 0) write(*,*) 'Could not allocate MnMinteraction!'
829    
830     end function MnMinitialize
831    
832     subroutine createInteractionLookup(this)
833     type(MnMinteractionMap), pointer :: this
834     integer :: biggestAtid, i, metal_atid, nonmetal_atid, error
835    
836     biggestAtid=-1
837     do i = 1, this%interactionCount
838     metal_atid = this%interactions(i)%metal_atid
839     nonmetal_atid = this%interactions(i)%nonmetal_atid
840    
841     if (metal_atid .gt. biggestAtid) biggestAtid = metal_atid
842     if (nonmetal_atid .gt. biggestAtid) biggestAtid = nonmetal_atid
843     enddo
844    
845     allocate(MnMinteractionLookup(biggestAtid,biggestAtid), stat=error)
846     if (error /= 0) write(*,*) 'Could not allocate MnMinteractionLookup'
847    
848     do i = 1, this%interactionCount
849     metal_atid = this%interactions(i)%metal_atid
850     nonmetal_atid = this%interactions(i)%nonmetal_atid
851    
852     MnMinteractionLookup(metal_atid, nonmetal_atid) = i
853     MnMinteractionLookup(nonmetal_atid, metal_atid) = i
854     enddo
855     end subroutine createInteractionLookup
856    
857    
858     function MnMdestroy(this) result(null_this)
859     logical :: done
860     type(MnMinteractionMap), pointer :: this
861     type(MnMinteractionMap), pointer :: null_this
862    
863     if (.not. associated(this)) then
864     null_this => null()
865     return
866     end if
867    
868     !! Walk down the list and deallocate each of the map's components
869     if(associated(this%interactions)) then
870     deallocate(this%interactions)
871     this%interactions=>null()
872     endif
873     deallocate(this)
874     this => null()
875     null_this => null()
876     end function MnMdestroy
877    
878    
879     subroutine deleteInteractions()
880     MnM_Map => MnMdestroy(MnM_Map)
881     return
882     end subroutine deleteInteractions
883    
884 chuckv 1173
885     subroutine getLJfunc(r, myPot, myDeriv)
886    
887     real(kind=dp), intent(in) :: r
888     real(kind=dp), intent(inout) :: myPot, myDeriv
889     real(kind=dp) :: ri, ri2, ri6, ri7, ri12, ri13
890     real(kind=dp) :: a, b, c, d, dx
891     integer :: j
892    
893     ri = 1.0_DP / r
894     ri2 = ri*ri
895     ri6 = ri2*ri2*ri2
896     ri7 = ri6*ri
897     ri12 = ri6*ri6
898     ri13 = ri12*ri
899    
900     myPot = 4.0_DP * (ri12 - ri6)
901     myDeriv = 24.0_DP * (ri7 - 2.0_DP * ri13)
902    
903     return
904     end subroutine getLJfunc
905    
906     subroutine getSoftFunc(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
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     myPot = 4.0_DP * (ri6)
919     myDeriv = - 24.0_DP * ri7
920    
921     return
922     end subroutine getSoftFunc
923    
924    
925    
926    
927    
928    
929 chuckv 1168 end module MetalNonMetal