ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 1231
Committed: Tue Mar 11 21:06:55 2008 UTC (17 years, 4 months ago) by chuckv
File size: 28538 byte(s)
Log Message:
MnM interaction now works.

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