ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/MetalNonMetal.F90
Revision: 3207
Committed: Tue Aug 14 17:41:05 2007 UTC (17 years, 11 months ago) by xsun
File size: 30124 byte(s)
Log Message:
added a check for non-initialized MnM_map

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.6 2007-08-14 17:41:05 xsun Exp $, $Date: 2007-08-14 17:41:05 $, $Name: not supported by cvs2svn $, $Revision: 1.6 $
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 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 public :: do_mnm_pair
109
110 contains
111
112
113 subroutine do_mnm_pair(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
114 Pot, A, F,t, Do_pot)
115 integer, intent(in) :: atom1, atom2
116 integer :: atid1, atid2, ljt1, ljt2
117 real( kind = dp ), intent(in) :: rij, r2, rcut
118 real( kind = dp ) :: pot, sw, vpair
119 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
120 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 real( kind = dp ), intent(inout), dimension(3) :: fpair
124 logical, intent(in) :: do_pot
125
126 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 if(.not.haveInteractionLookup) then
138 call createInteractionLookup(MnM_MAP)
139 haveInteractionLookup =.true.
140 end if
141
142 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 call calc_mnm_lennardjones(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, &
148 Fpair, Pot, F, Do_pot, interaction_id)
149 case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE)
150 call calc_mnm_morse(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
151 Pot, F, Do_pot, interaction_id, interaction_type)
152 case(MNM_MAW)
153 call calc_mnm_maw(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, &
154 Pot,A, F,t, Do_pot, interaction_id)
155 end select
156
157 end subroutine do_mnm_pair
158
159 subroutine calc_mnm_lennardjones(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, &
160 Fpair, Pot, F, Do_pot, interaction_id)
161
162 integer, intent(in) :: atom1, atom2
163 real( kind = dp ), intent(in) :: rij, r2, rcut
164 real( kind = dp ) :: pot, sw, vpair
165 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
166 real( kind = dp ), intent(in), dimension(3) :: d
167 real( kind = dp ), intent(inout), dimension(3) :: fpair
168 logical, intent(in) :: do_pot
169 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 sigmai = 1.0_dp / MnM_Map%interactions(interaction_id)%sigma
182 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 endif
202
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
215 #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
221 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
225 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
229 #else
230 if (do_pot) pot = pot + sw*pot_temp
231
232 f(1,atom1) = f(1,atom1) + fx
233 f(2,atom1) = f(2,atom1) + fy
234 f(3,atom1) = f(3,atom1) + fz
235
236 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
249 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
250
251 fpair(1) = fpair(1) + fx
252 fpair(2) = fpair(2) + fy
253 fpair(3) = fpair(3) + fz
254
255 endif
256 return
257 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 integer, intent(in) :: atom1, atom2
262 real( kind = dp ), intent(in) :: rij, r2, rcut
263 real( kind = dp ) :: pot, sw, vpair
264 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
265 real( kind = dp ), intent(in), dimension(3) :: d
266 real( kind = dp ), intent(inout), dimension(3) :: fpair
267 logical, intent(in) :: do_pot
268 integer, intent(in) :: interaction_id, interaction_type
269 logical :: shiftedPot, shiftedFrc
270
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
299
300 D0 = MnM_Map%interactions(interaction_id)%D0
301 R0 = MnM_Map%interactions(interaction_id)%r0
302 Beta0 = MnM_Map%interactions(interaction_id)%Beta0
303 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 expfnc = exp(expt)
310 expfnc2 = expfnc*expfnc
311
312 if (shiftedPot .or. shiftedFrc) then
313 exptC = -beta0*(rcut - R0)
314 expfncC = exp(exptC)
315 expfnc2C = expfncC*expfncC
316 endif
317
318 select case (interaction_type)
319 case (MNM_SHIFTEDMORSE)
320
321 myPot = D0 * (expfnc2 - 2.0_dp * expfnc)
322 myDeriv = 2.0*D0*beta0*(expfnc - expfnc2)
323
324 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
336 case (MNM_REPULSIVEMORSE)
337
338 myPot = D0 * expfnc2
339 myDeriv = -2.0_dp * D0 * beta0 * expfnc2
340
341 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 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 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 real( kind = dp ), intent(inout), dimension(3,nLocal) :: f
417 real (kind=dp),intent(in), dimension(9,nLocal) :: A
418 real (kind=dp),intent(inout), dimension(3,nLocal) :: t
419
420 real( kind = dp ), intent(in), dimension(3) :: d
421 real( kind = dp ), intent(inout), dimension(3) :: fpair
422 logical, intent(in) :: do_pot
423
424 integer, intent(in) :: interaction_id
425
426 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
440
441
442 #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
465 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
483 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
490
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
503
504 !!$ 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 r3 = r2*rij
524 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 ! These derivatives can be obtained by using
530 ! cos(theta) = \hat{u} . \vec{r} / r
531 ! where \hat{u} is the body-fixed unit frame for the water molecule,
532 ! and \vec{r} is the vector to the metal atom.
533 !
534 ! the derivatives wrt \vec{r} are:
535 ! dcostheta/d\vec{r} = - costheta \vec{r} / r^2 + \hat{u} / r
536 ! the molecular frame for each water has u = {0, 0, 1}, so these:
537 !
538 ! dctdx = - x * z / r3 + ux / rij
539 ! dctdy = - y * z / r3 + uy / rij
540 ! dctdz = - z2 / r3 + uz / rij
541 !
542 ! become:
543 !
544 dctdx = - z * x / r3
545 dctdy = - z * y / r3
546 dctdz = 1.0_dp / rij - z2 / r3
547
548 dctdux = x / rij
549 dctduy = y / rij
550 dctduz = z / rij
551
552 ! dctdux = y / rij ! - (z * x2) / r3
553 ! dctduy = -x / rij !- (z * y2) / r3
554 ! dctduz = 0.0_dp ! z / rij - (z2 * z) / r3
555
556 ! this is an attempt to try to truncate the singularity when
557 ! sin(theta) is near 0.0:
558
559 st2 = 1.0_dp - ct*ct
560 if (abs(st2) .lt. 1.0e-12_dp) then
561 proj = sqrt(rij * 1.0e-12_dp)
562 !!$ dcpdx = 1.0_dp / proj
563 !!$ dcpdy = 0.0_dp
564 !!$ dcpdux = x / proj
565 !!$ dcpduy = 0.0_dp
566 dspdx = 0.0_dp
567 dspdy = 1.0_dp / proj
568 dspdux = 0.0_dp
569 dspduy = y / proj
570 else
571 proj = sqrt(x2 + y2)
572 proj3 = proj*proj*proj
573 !!$ dcpdx = 1.0_dp / proj - x2 / proj3
574 !!$ dcpdy = - x * y / proj3
575 !!$ dcpdux = x / proj - (x2 * x) / proj3
576 !!$ dcpduy = - (x * y2) / proj3
577 dspdx = - x * y / proj3
578 dspdy = 1.0_dp / proj - y2 / proj3
579 dspdux = - (y * x2) / proj3
580 dspduy = y / proj - (y2 * y) / proj3
581 endif
582
583 !!$ cp = x / proj
584 !!$ dcpdz = 0.0_dp
585 !!$ dcpduz = 0.0_dp
586
587 sp = y / proj
588 dspdz = 0.0_dp
589 dspduz = 0.0_dp
590
591
592 pot_temp = D0 * (expfnc02 - 2.0_dp * expfnc0) + &
593 2.0_dp * gamma * D0 * expfncH2 * (1.0_dp + alpha*ct)*sp*sp
594 vpair = vpair + pot_temp
595
596 if (do_pot) then
597 #ifdef IS_MPI
598 pot_row(VDW_POT,atom1) = pot_row(VDW_POT,atom1) + 0.5_dp*pot_temp*sw
599 pot_col(VDW_POT,atom2) = pot_col(VDW_POT,atom2) + 0.5_dp*pot_temp*sw
600 #else
601 pot = pot + pot_temp*sw
602 #endif
603 endif
604
605 ! derivative wrt r
606 dmawdr = 2.0*D0*beta0*(expfnc0 - expfnc02) - &
607 4.0 * gamma * D0 * betaH * expfncH2 * (1.0_dp + alpha*ct)*sp*sp
608
609 ! derivative wrt cos(theta)
610 dmawdct = 2.0 * gamma * D0 * expfncH2 * alpha * sp * sp
611
612 ! derivative wrt sin(phi)
613 dmawdsp = 4.0 * gamma * D0 * expfncH2 * (1.0_dp + alpha*ct)*sp
614
615
616 ! chain rule to put these back on x, y, z, ux, uy, uz
617 dvdx = dmawdr * drdx + dmawdct * dctdx + dmawdsp * dspdx
618 dvdy = dmawdr * drdy + dmawdct * dctdy + dmawdsp * dspdy
619 dvdz = dmawdr * drdz + dmawdct * dctdz + dmawdsp * dspdz
620
621
622
623
624 dvdux = dmawdr * drdux + dmawdct * dctdux + dmawdsp * dspdux
625 dvduy = dmawdr * drduy + dmawdct * dctduy + dmawdsp * dspduy
626 dvduz = dmawdr * drduz + dmawdct * dctduz + dmawdsp * dspduz
627
628
629
630
631 tx = (dvduy - dvduz) * sw
632 ty = (dvduz - dvdux) * sw
633 tz = (dvdux - dvduy) * sw
634
635 ! go back to lab frame using transpose of rotation matrix:
636
637 #ifdef IS_MPI
638 if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
639 t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*tx + &
640 a_Row(4,atom1)*ty + a_Row(7,atom1)*tz
641 t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*tx + &
642 a_Row(5,atom1)*ty + a_Row(8,atom1)*tz
643 t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*tx + &
644 a_Row(6,atom1)*ty + a_Row(9,atom1)*tz
645 else
646 t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*tx + &
647 a_Col(4,atom2)*ty + a_Col(7,atom2)*tz
648 t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*tx + &
649 a_Col(5,atom2)*ty + a_Col(8,atom2)*tz
650 t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*tx + &
651 a_Col(6,atom2)*ty + a_Col(9,atom2)*tz
652 endif
653 #else
654 if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
655 t(1,atom1) = t(1,atom1) + a(1,atom1)*tx + a(4,atom1)*ty + &
656 a(7,atom1)*tz
657 t(2,atom1) = t(2,atom1) + a(2,atom1)*tx + a(5,atom1)*ty + &
658 a(8,atom1)*tz
659 t(3,atom1) = t(3,atom1) + a(3,atom1)*tx + a(6,atom1)*ty + &
660 a(9,atom1)*tz
661 else
662 t(1,atom2) = t(1,atom2) + a(1,atom2)*tx + a(4,atom2)*ty + &
663 a(7,atom2)*tz
664 t(2,atom2) = t(2,atom2) + a(2,atom2)*tx + a(5,atom2)*ty + &
665 a(8,atom2)*tz
666 t(3,atom2) = t(3,atom2) + a(3,atom2)*tx + a(6,atom2)*ty + &
667 a(9,atom2)*tz
668 endif
669 #endif
670 ! Now, on to the forces:
671
672 fx = -dvdx * sw
673 fy = -dvdy * sw
674 fz = -dvdz * sw
675
676 ! rotate the terms back into the lab frame:
677
678 #ifdef IS_MPI
679 if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
680 fxl = a_Row(1,atom1)*fx + a_Row(4,atom1)*fy + a_Row(7,atom1)*fz
681 fyl = a_Row(2,atom1)*fx + a_Row(5,atom1)*fy + a_Row(8,atom1)*fz
682 fzl = a_Row(3,atom1)*fx + a_Row(6,atom1)*fy + a_Row(9,atom1)*fz
683 else
684 fxl = a_Col(1,atom2)*fx + a_Col(4,atom2)*fy + a_Col(7,atom2)*fz
685 fyl = a_Col(2,atom2)*fx + a_Col(5,atom2)*fy + a_Col(8,atom2)*fz
686 fzl = a_Col(3,atom2)*fx + a_Col(6,atom2)*fy + a_Col(9,atom2)*fz
687 endif
688 f_Row(1,atom1) = f_Row(1,atom1) + fxl
689 f_Row(2,atom1) = f_Row(2,atom1) + fyl
690 f_Row(3,atom1) = f_Row(3,atom1) + fzl
691
692 f_Col(1,atom2) = f_Col(1,atom2) - fxl
693 f_Col(2,atom2) = f_Col(2,atom2) - fyl
694 f_Col(3,atom2) = f_Col(3,atom2) - fzl
695 #else
696 if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then
697 fxl = a(1,atom1)*fx + a(4,atom1)*fy + a(7,atom1)*fz
698 fyl = a(2,atom1)*fx + a(5,atom1)*fy + a(8,atom1)*fz
699 fzl = a(3,atom1)*fx + a(6,atom1)*fy + a(9,atom1)*fz
700 else
701 fxl = a(1,atom2)*fx + a(4,atom2)*fy + a(7,atom2)*fz
702 fyl = a(2,atom2)*fx + a(5,atom2)*fy + a(8,atom2)*fz
703 fzl = a(3,atom2)*fx + a(6,atom2)*fy + a(9,atom2)*fz
704 endif
705 f(1,atom1) = f(1,atom1) + fxl
706 f(2,atom1) = f(2,atom1) + fyl
707 f(3,atom1) = f(3,atom1) + fzl
708
709 f(1,atom2) = f(1,atom2) - fxl
710 f(2,atom2) = f(2,atom2) - fyl
711 f(3,atom2) = f(3,atom2) - fzl
712 #endif
713
714 #ifdef IS_MPI
715 id1 = AtomRowToGlobal(atom1)
716 id2 = AtomColToGlobal(atom2)
717 #else
718 id1 = atom1
719 id2 = atom2
720 #endif
721
722 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
723
724 fpair(1) = fpair(1) + fxl
725 fpair(2) = fpair(2) + fyl
726 fpair(3) = fpair(3) + fzl
727
728 endif
729
730 return
731 end subroutine calc_mnm_maw
732
733
734 subroutine setMnMDefaultCutoff(thisRcut, shiftedPot, shiftedFrc)
735 real(kind=dp), intent(in) :: thisRcut
736 logical, intent(in) :: shiftedPot
737 logical, intent(in) :: shiftedFrc
738 integer i, nInteractions
739 defaultCutoff = thisRcut
740 defaultShiftPot = shiftedPot
741 defaultShiftFrc = shiftedFrc
742
743 if (associated(MnM_Map)) then
744 if(MnM_Map%interactionCount /= 0) then
745 nInteractions = MnM_Map%interactionCount
746
747 do i = 1, nInteractions
748 MnM_Map%interactions(i)%shiftedPot = shiftedPot
749 MnM_Map%interactions(i)%shiftedFrc = shiftedFrc
750 MnM_Map%interactions(i)%rCut = thisRcut
751 MnM_Map%interactions(i)%rCutWasSet = .true.
752 enddo
753 end if
754 end if
755
756 end subroutine setMnMDefaultCutoff
757
758 subroutine copyAllData(v1, v2)
759 type(MnMinteractionMap), pointer :: v1
760 type(MnMinteractionMap), pointer :: v2
761 integer :: i, j
762
763 do i = 1, v1%interactionCount
764 v2%interactions(i) = v1%interactions(i)
765 enddo
766
767 v2%interactionCount = v1%interactionCount
768 return
769 end subroutine copyAllData
770
771 subroutine addInteraction(myInteraction)
772 type(MNMtype) :: myInteraction
773 type(MnMinteraction) :: nt
774 integer :: id
775
776
777
778
779
780 nt%interaction_type = myInteraction%MNMInteractionType
781 nt%metal_atid = &
782 getFirstMatchingElement(atypes, "c_ident", myInteraction%metal_atid)
783 nt%nonmetal_atid = &
784 getFirstMatchingElement(atypes, "c_ident", myInteraction%nonmetal_atid)
785
786
787 select case (nt%interaction_type)
788 case (MNM_LENNARDJONES)
789 nt%sigma = myInteraction%sigma
790 nt%epsilon = myInteraction%epsilon
791 case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE)
792 nt%R0 = myInteraction%R0
793 nt%D0 = myInteraction%D0
794 nt%beta0 = myInteraction%beta0
795 case(MNM_MAW)
796 nt%R0 = myInteraction%R0
797 nt%D0 = myInteraction%D0
798 nt%beta0 = myInteraction%beta0
799 nt%betaH = myInteraction%betaH
800 nt%alpha = myInteraction%alpha
801 nt%gamma = myInteraction%gamma
802 case default
803 call handleError("MNM", "Unknown Interaction type")
804 end select
805
806 if (.not. associated(MnM_Map)) then
807 call ensureCapacityHelper(MnM_Map, 1)
808 else
809 call ensureCapacityHelper(MnM_Map, MnM_Map%interactionCount + 1)
810 end if
811
812 MnM_Map%interactionCount = MnM_Map%interactionCount + 1
813 id = MnM_Map%interactionCount
814 MnM_Map%interactions(id) = nt
815 end subroutine addInteraction
816
817 subroutine ensureCapacityHelper(this, minCapacity)
818 type(MnMinteractionMap), pointer :: this, that
819 integer, intent(in) :: minCapacity
820 integer :: oldCapacity
821 integer :: newCapacity
822 logical :: resizeFlag
823
824 resizeFlag = .false.
825
826 ! first time: allocate a new vector with default size
827
828 if (.not. associated(this)) then
829 this => MnMinitialize(minCapacity, 0)
830 endif
831
832 oldCapacity = size(this%interactions)
833
834 if (minCapacity > oldCapacity) then
835 if (this%capacityIncrement .gt. 0) then
836 newCapacity = oldCapacity + this%capacityIncrement
837 else
838 newCapacity = oldCapacity * 2
839 endif
840 if (newCapacity .lt. minCapacity) then
841 newCapacity = minCapacity
842 endif
843 resizeFlag = .true.
844 else
845 newCapacity = oldCapacity
846 endif
847
848 if (resizeFlag) then
849 that => MnMinitialize(newCapacity, this%capacityIncrement)
850 call copyAllData(this, that)
851 this => MnMdestroy(this)
852 this => that
853 endif
854 end subroutine ensureCapacityHelper
855
856 function MnMinitialize(cap, capinc) result(this)
857 integer, intent(in) :: cap, capinc
858 integer :: error
859 type(MnMinteractionMap), pointer :: this
860
861 nullify(this)
862
863 if (cap < 0) then
864 write(*,*) 'Bogus Capacity:', cap
865 return
866 endif
867 allocate(this,stat=error)
868 if ( error /= 0 ) then
869 write(*,*) 'Could not allocate MnMinteractionMap!'
870 return
871 end if
872
873 this%initialCapacity = cap
874 this%capacityIncrement = capinc
875
876 allocate(this%interactions(this%initialCapacity), stat=error)
877 if(error /= 0) write(*,*) 'Could not allocate MnMinteraction!'
878
879 end function MnMinitialize
880
881 subroutine createInteractionLookup(this)
882 type (MnMInteractionMap),pointer :: this
883 integer :: biggestAtid, i, metal_atid, nonmetal_atid, error
884
885 biggestAtid=-1
886 do i = 1, this%interactionCount
887 metal_atid = this%interactions(i)%metal_atid
888 nonmetal_atid = this%interactions(i)%nonmetal_atid
889
890 if (metal_atid .gt. biggestAtid) biggestAtid = metal_atid
891 if (nonmetal_atid .gt. biggestAtid) biggestAtid = nonmetal_atid
892 enddo
893
894 allocate(MnMinteractionLookup(biggestAtid,biggestAtid), stat=error)
895 if (error /= 0) write(*,*) 'Could not allocate MnMinteractionLookup'
896
897 do i = 1, this%interactionCount
898 metal_atid = this%interactions(i)%metal_atid
899 nonmetal_atid = this%interactions(i)%nonmetal_atid
900
901 MnMinteractionLookup(metal_atid, nonmetal_atid) = i
902 MnMinteractionLookup(nonmetal_atid, metal_atid) = i
903 enddo
904 end subroutine createInteractionLookup
905
906
907 function MnMdestroy(this) result(null_this)
908 logical :: done
909 type(MnMinteractionMap), pointer :: this
910 type(MnMinteractionMap), pointer :: null_this
911
912 if (.not. associated(this)) then
913 null_this => null()
914 return
915 end if
916
917 !! Walk down the list and deallocate each of the map's components
918 if(associated(this%interactions)) then
919 deallocate(this%interactions)
920 this%interactions=>null()
921 endif
922 deallocate(this)
923 this => null()
924 null_this => null()
925 end function MnMdestroy
926
927
928 subroutine deleteInteractions()
929 MnM_Map => MnMdestroy(MnM_Map)
930 return
931 end subroutine deleteInteractions
932
933
934 subroutine getLJfunc(r, myPot, myDeriv)
935
936 real(kind=dp), intent(in) :: r
937 real(kind=dp), intent(inout) :: myPot, myDeriv
938 real(kind=dp) :: ri, ri2, ri6, ri7, ri12, ri13
939 real(kind=dp) :: a, b, c, d, dx
940 integer :: j
941
942 ri = 1.0_DP / r
943 ri2 = ri*ri
944 ri6 = ri2*ri2*ri2
945 ri7 = ri6*ri
946 ri12 = ri6*ri6
947 ri13 = ri12*ri
948
949 myPot = 4.0_DP * (ri12 - ri6)
950 myDeriv = 24.0_DP * (ri7 - 2.0_DP * ri13)
951
952 return
953 end subroutine getLJfunc
954
955 subroutine getSoftFunc(r, myPot, myDeriv)
956
957 real(kind=dp), intent(in) :: r
958 real(kind=dp), intent(inout) :: myPot, myDeriv
959 real(kind=dp) :: ri, ri2, ri6, ri7
960 real(kind=dp) :: a, b, c, d, dx
961 integer :: j
962
963 ri = 1.0_DP / r
964 ri2 = ri*ri
965 ri6 = ri2*ri2*ri2
966 ri7 = ri6*ri
967 myPot = 4.0_DP * (ri6)
968 myDeriv = - 24.0_DP * ri7
969
970 return
971 end subroutine getSoftFunc
972
973
974
975
976
977
978 end module MetalNonMetal