ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 1175
Committed: Thu Jul 19 19:49:38 2007 UTC (18 years ago) by chuckv
Original Path: trunk/src/UseTheForce/DarkSide/MetalNonMetal.F90
File size: 29450 byte(s)
Log Message:
More bug fixes to Metalnonmetal.

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