ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 1388
Committed: Fri Oct 30 16:38:48 2009 UTC (15 years, 9 months ago) by chuckv
File size: 24727 byte(s)
Log Message:
Removed remaining MPI from metallic potentials. Bug in MPI calculation of energies in Sutton-Chen.

File Contents

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