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.8 2008-03-03 17:14:36 chuckv Exp $, $Date: 2008-03-03 17:14:36 $, $Name: not supported by cvs2svn $, $Revision: 1.8 $ |
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, proj, proj3, st2,st,s2t |
435 |
real(kind = dp) :: drdx, drdy, drdz, drdux, drduy, drduz |
436 |
real(kind = dp) :: ct,ct2,c2t, dctdx, dctdy, dctdz, dctdux, dctduy, dctduz |
437 |
real(kind = dp) :: sp,c2p,sp2,dspdx, dspdy, dspdz, dspdux, dspduy, dspduz |
438 |
real(kind = dp) :: dvdx, dvdy, dvdz, dvdux, dvduy, dvduz |
439 |
real(kind = dp) :: Vang, dVangdx, dVangdy, dVangdz, dVangdux, dVangduy, dVangduz |
440 |
real(kind = dp) :: dVangdct, dVangdsp, dVmorsedr |
441 |
real(kind = dp) :: Vmorse, dVmorsedx, dVmorsedy, dVmorsedz |
442 |
real(kind = dp) :: dVmorsedux, dVmorseduy, dVmorseduz |
443 |
real(kind = dp) :: fx, fy, fz, tx, ty, tz, fxl, fyl, fzl |
444 |
integer :: atid1, atid2, id1, id2 |
445 |
logical :: shiftedPot, shiftedFrc |
446 |
|
447 |
|
448 |
|
449 |
|
450 |
#ifdef IS_MPI |
451 |
atid1 = atid_Row(atom1) |
452 |
atid2 = atid_Col(atom2) |
453 |
|
454 |
if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then |
455 |
! rotate the inter-particle separation into the two different |
456 |
! body-fixed coordinate systems: |
457 |
|
458 |
x = A_row(1,atom1)*d(1) + A_row(2,atom1)*d(2) + A_row(3,atom1)*d(3) |
459 |
y = A_row(4,atom1)*d(1) + A_row(5,atom1)*d(2) + A_row(6,atom1)*d(3) |
460 |
z = A_row(7,atom1)*d(1) + A_row(8,atom1)*d(2) + A_row(9,atom1)*d(3) |
461 |
else |
462 |
! negative sign because this is the vector from j to i: |
463 |
|
464 |
x = -(A_Col(1,atom2)*d(1) + A_Col(2,atom2)*d(2) + A_Col(3,atom2)*d(3)) |
465 |
y = -(A_Col(4,atom2)*d(1) + A_Col(5,atom2)*d(2) + A_Col(6,atom2)*d(3)) |
466 |
z = -(A_Col(7,atom2)*d(1) + A_Col(8,atom2)*d(2) + A_Col(9,atom2)*d(3)) |
467 |
endif |
468 |
|
469 |
#else |
470 |
atid1 = atid(atom1) |
471 |
atid2 = atid(atom2) |
472 |
|
473 |
if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then |
474 |
! rotate the inter-particle separation into the two different |
475 |
! body-fixed coordinate systems: |
476 |
|
477 |
x = a(1,atom1)*d(1) + a(2,atom1)*d(2) + a(3,atom1)*d(3) |
478 |
y = a(4,atom1)*d(1) + a(5,atom1)*d(2) + a(6,atom1)*d(3) |
479 |
z = a(7,atom1)*d(1) + a(8,atom1)*d(2) + a(9,atom1)*d(3) |
480 |
else |
481 |
! negative sign because this is the vector from j to i: |
482 |
|
483 |
x = -(a(1,atom2)*d(1) + a(2,atom2)*d(2) + a(3,atom2)*d(3)) |
484 |
y = -(a(4,atom2)*d(1) + a(5,atom2)*d(2) + a(6,atom2)*d(3)) |
485 |
z = -(a(7,atom2)*d(1) + a(8,atom2)*d(2) + a(9,atom2)*d(3)) |
486 |
endif |
487 |
|
488 |
#endif |
489 |
|
490 |
|
491 |
D0 = MnM_Map%interactions(interaction_id)%D0 |
492 |
R0 = MnM_Map%interactions(interaction_id)%r0 |
493 |
beta0 = MnM_Map%interactions(interaction_id)%beta0 |
494 |
alpha = MnM_Map%interactions(interaction_id)%alpha |
495 |
|
496 |
|
497 |
|
498 |
|
499 |
shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot |
500 |
shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc |
501 |
|
502 |
expt0 = -beta0*(rij - R0) |
503 |
expfnc0 = exp(expt0) |
504 |
expfnc02 = expfnc0*expfnc0 |
505 |
|
506 |
|
507 |
|
508 |
|
509 |
!!$ if (shiftedPot .or. shiftedFrc) then |
510 |
!!$ exptC0 = -beta0*(rcut - R0) |
511 |
!!$ expfncC0 = exp(exptC0) |
512 |
!!$ expfncC02 = expfncC0*expfncC0 |
513 |
!!$ exptCH = -betaH*(rcut - R0) |
514 |
!!$ expfncCH = exp(exptCH) |
515 |
!!$ expfncCH2 = expfncCH*expfncCH |
516 |
!!$ endif |
517 |
|
518 |
drdx = -d(1) / rij |
519 |
drdy = -d(2) / rij |
520 |
drdz = -d(3) / rij |
521 |
drdux = 0.0_dp |
522 |
drduy = 0.0_dp |
523 |
drduz = 0.0_dp |
524 |
|
525 |
x2 = x*x |
526 |
y2 = y*y |
527 |
z2 = z*z |
528 |
r3 = r2*rij |
529 |
ct = z / rij |
530 |
ct2 = ct * ct |
531 |
|
532 |
if (ct .gt. 1.0_dp) ct = 1.0_dp |
533 |
if (ct .lt. -1.0_dp) ct = -1.0_dp |
534 |
|
535 |
! These derivatives can be obtained by using |
536 |
! cos(theta) = \hat{u} . \vec{r} / r |
537 |
! where \hat{u} is the body-fixed unit frame for the water molecule, |
538 |
! and \vec{r} is the vector to the metal atom. |
539 |
! |
540 |
! the derivatives wrt \vec{r} are: |
541 |
! dcostheta/d\vec{r} = - costheta \vec{r} / r^2 + \hat{u} / r |
542 |
! the molecular frame for each water has u = {0, 0, 1}, so these: |
543 |
! |
544 |
! dctdx = - x * z / r3 + ux / rij |
545 |
! dctdy = - y * z / r3 + uy / rij |
546 |
! dctdz = - z2 / r3 + uz / rij |
547 |
! |
548 |
! become: |
549 |
! |
550 |
dctdx = - z * x / r3 |
551 |
dctdy = - z * y / r3 |
552 |
dctdz = 1.0_dp / rij - z2 / r3 |
553 |
|
554 |
dctdux = x / rij |
555 |
dctduy = y / rij |
556 |
dctduz = z / rij |
557 |
|
558 |
! this is an attempt to try to truncate the singularity when |
559 |
! sin(theta) is near 0.0: |
560 |
|
561 |
st2 = 1.0_dp - ct2 |
562 |
|
563 |
if (abs(st2) .lt. 1.0e-12_dp) then |
564 |
proj = sqrt(rij * 1.0e-12_dp) |
565 |
dspdx = 0.0_dp |
566 |
dspdy = 1.0_dp / proj |
567 |
dspdux = 0.0_dp |
568 |
dspduy = y / proj |
569 |
else |
570 |
proj = sqrt(x2 + y2) |
571 |
proj3 = proj*proj*proj |
572 |
dspdx = - x * y / proj3 |
573 |
dspdy = 1.0_dp / proj - y2 / proj3 |
574 |
dspdux = - (y * x2) / proj3 |
575 |
dspduy = y / proj - (y2 * y) / proj3 |
576 |
endif |
577 |
|
578 |
sp = y / proj |
579 |
dspdz = 0.0_dp |
580 |
dspduz = 0.0_dp |
581 |
|
582 |
sp2 = sp * sp |
583 |
c2p = 1.0_dp - 2.0_dp * sp2 |
584 |
|
585 |
! V(r) = D_e exp(-a(r-re)(exp(-a(r-re))-2) |
586 |
Vmorse = D0 * (expfnc02 - 2.0_dp * expfnc0) |
587 |
! angular modulation of morse part of potential to approximate a sp3 orbital |
588 |
! Vang = 1 - alpha*(1/2 + sqrt(3.0)*cos(theta)/4.0 - 3.0*cos(2.0*phi)*sin^2(theta)/8.0) |
589 |
|
590 |
Vang = 1.0_dp - alpha*(0.5_dp + sqrt(3.0_dp)*ct/4.0_dp - 3.0*c2p*st2/8.0_dp) |
591 |
|
592 |
pot_temp = Vmorse*Vang |
593 |
|
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 |
dVmorsedr = 2.0_dp*D0*beta0*(expfnc0 - expfnc02) |
606 |
dVmorsedx = dVmorsedr * drdx |
607 |
dVmorsedy = dVmorsedr * drdy |
608 |
dVmorsedz = dVmorsedr * drdz |
609 |
dVmorsedux = 0.0_dp |
610 |
dVmorseduy = 0.0_dp |
611 |
dVmorseduy = 0.0_dp |
612 |
|
613 |
dVangdct = -alpha * ( sqrt(3.0_dp) + 3*ct*c2p ) / 4.0_dp |
614 |
dVangdsp = -3.0_dp * alpha * st2 * sp / 2.0_dp |
615 |
|
616 |
dVangdx = dVangdct * dctdx + dVangdsp * dspdx |
617 |
dVangdy = dVangdct * dctdy + dVangdsp * dspdy |
618 |
dVangdy = dVangdct * dctdy + dVangdsp * dspdy |
619 |
dVangdux = dVangdct * dctdux + dVangdsp * dspdux |
620 |
dVangduy = dVangdct * dctduy + dVangdsp * dspduy |
621 |
dVangduy = dVangdct * dctduy + dVangdsp * dspduy |
622 |
|
623 |
! chain rule to put these back on x, y, z, ux, uy, uz |
624 |
dvdx = Vang * dVmorsedx + Vmorse * dVangdx |
625 |
dvdy = Vang * dVmorsedy + Vmorse * dVangdy |
626 |
dvdz = Vang * dVmorsedz + Vmorse * dVangdz |
627 |
|
628 |
dvdux = Vang * dVmorsedux + Vmorse * dVangdux |
629 |
dvduy = Vang * dVmorseduy + Vmorse * dVangduy |
630 |
dvduz = Vang * dVmorseduz + Vmorse * dVangduz |
631 |
|
632 |
tx = (dvduy - dvduz) * sw |
633 |
ty = (dvduz - dvdux) * sw |
634 |
tz = (dvdux - dvduy) * sw |
635 |
|
636 |
! go back to lab frame using transpose of rotation matrix: |
637 |
|
638 |
#ifdef IS_MPI |
639 |
if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then |
640 |
t_Row(1,atom1) = t_Row(1,atom1) + a_Row(1,atom1)*tx + & |
641 |
a_Row(4,atom1)*ty + a_Row(7,atom1)*tz |
642 |
t_Row(2,atom1) = t_Row(2,atom1) + a_Row(2,atom1)*tx + & |
643 |
a_Row(5,atom1)*ty + a_Row(8,atom1)*tz |
644 |
t_Row(3,atom1) = t_Row(3,atom1) + a_Row(3,atom1)*tx + & |
645 |
a_Row(6,atom1)*ty + a_Row(9,atom1)*tz |
646 |
else |
647 |
t_Col(1,atom2) = t_Col(1,atom2) + a_Col(1,atom2)*tx + & |
648 |
a_Col(4,atom2)*ty + a_Col(7,atom2)*tz |
649 |
t_Col(2,atom2) = t_Col(2,atom2) + a_Col(2,atom2)*tx + & |
650 |
a_Col(5,atom2)*ty + a_Col(8,atom2)*tz |
651 |
t_Col(3,atom2) = t_Col(3,atom2) + a_Col(3,atom2)*tx + & |
652 |
a_Col(6,atom2)*ty + a_Col(9,atom2)*tz |
653 |
endif |
654 |
#else |
655 |
if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then |
656 |
t(1,atom1) = t(1,atom1) + a(1,atom1)*tx + a(4,atom1)*ty + & |
657 |
a(7,atom1)*tz |
658 |
t(2,atom1) = t(2,atom1) + a(2,atom1)*tx + a(5,atom1)*ty + & |
659 |
a(8,atom1)*tz |
660 |
t(3,atom1) = t(3,atom1) + a(3,atom1)*tx + a(6,atom1)*ty + & |
661 |
a(9,atom1)*tz |
662 |
else |
663 |
t(1,atom2) = t(1,atom2) + a(1,atom2)*tx + a(4,atom2)*ty + & |
664 |
a(7,atom2)*tz |
665 |
t(2,atom2) = t(2,atom2) + a(2,atom2)*tx + a(5,atom2)*ty + & |
666 |
a(8,atom2)*tz |
667 |
t(3,atom2) = t(3,atom2) + a(3,atom2)*tx + a(6,atom2)*ty + & |
668 |
a(9,atom2)*tz |
669 |
endif |
670 |
#endif |
671 |
! Now, on to the forces: |
672 |
|
673 |
fx = dvdx * sw |
674 |
fy = dvdy * sw |
675 |
fz = dvdz * sw |
676 |
|
677 |
! rotate the terms back into the lab frame: |
678 |
|
679 |
#ifdef IS_MPI |
680 |
if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then |
681 |
fxl = a_Row(1,atom1)*fx + a_Row(4,atom1)*fy + a_Row(7,atom1)*fz |
682 |
fyl = a_Row(2,atom1)*fx + a_Row(5,atom1)*fy + a_Row(8,atom1)*fz |
683 |
fzl = a_Row(3,atom1)*fx + a_Row(6,atom1)*fy + a_Row(9,atom1)*fz |
684 |
else |
685 |
fxl = a_Col(1,atom2)*fx + a_Col(4,atom2)*fy + a_Col(7,atom2)*fz |
686 |
fyl = a_Col(2,atom2)*fx + a_Col(5,atom2)*fy + a_Col(8,atom2)*fz |
687 |
fzl = a_Col(3,atom2)*fx + a_Col(6,atom2)*fy + a_Col(9,atom2)*fz |
688 |
endif |
689 |
f_Row(1,atom1) = f_Row(1,atom1) + fxl |
690 |
f_Row(2,atom1) = f_Row(2,atom1) + fyl |
691 |
f_Row(3,atom1) = f_Row(3,atom1) + fzl |
692 |
|
693 |
f_Col(1,atom2) = f_Col(1,atom2) - fxl |
694 |
f_Col(2,atom2) = f_Col(2,atom2) - fyl |
695 |
f_Col(3,atom2) = f_Col(3,atom2) - fzl |
696 |
#else |
697 |
if (atid2.eq.MnM_Map%interactions(interaction_id)%metal_atid) then |
698 |
fxl = a(1,atom1)*fx + a(4,atom1)*fy + a(7,atom1)*fz |
699 |
fyl = a(2,atom1)*fx + a(5,atom1)*fy + a(8,atom1)*fz |
700 |
fzl = a(3,atom1)*fx + a(6,atom1)*fy + a(9,atom1)*fz |
701 |
else |
702 |
fxl = a(1,atom2)*fx + a(4,atom2)*fy + a(7,atom2)*fz |
703 |
fyl = a(2,atom2)*fx + a(5,atom2)*fy + a(8,atom2)*fz |
704 |
fzl = a(3,atom2)*fx + a(6,atom2)*fy + a(9,atom2)*fz |
705 |
endif |
706 |
f(1,atom1) = f(1,atom1) + fxl |
707 |
f(2,atom1) = f(2,atom1) + fyl |
708 |
f(3,atom1) = f(3,atom1) + fzl |
709 |
|
710 |
f(1,atom2) = f(1,atom2) - fxl |
711 |
f(2,atom2) = f(2,atom2) - fyl |
712 |
f(3,atom2) = f(3,atom2) - fzl |
713 |
#endif |
714 |
|
715 |
#ifdef IS_MPI |
716 |
id1 = AtomRowToGlobal(atom1) |
717 |
id2 = AtomColToGlobal(atom2) |
718 |
#else |
719 |
id1 = atom1 |
720 |
id2 = atom2 |
721 |
#endif |
722 |
|
723 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
724 |
|
725 |
fpair(1) = fpair(1) + fxl |
726 |
fpair(2) = fpair(2) + fyl |
727 |
fpair(3) = fpair(3) + fzl |
728 |
|
729 |
endif |
730 |
|
731 |
return |
732 |
end subroutine calc_mnm_maw |
733 |
|
734 |
|
735 |
subroutine setMnMDefaultCutoff(thisRcut, shiftedPot, shiftedFrc) |
736 |
real(kind=dp), intent(in) :: thisRcut |
737 |
logical, intent(in) :: shiftedPot |
738 |
logical, intent(in) :: shiftedFrc |
739 |
integer i, nInteractions |
740 |
defaultCutoff = thisRcut |
741 |
defaultShiftPot = shiftedPot |
742 |
defaultShiftFrc = shiftedFrc |
743 |
|
744 |
if (associated(MnM_Map)) then |
745 |
if(MnM_Map%interactionCount /= 0) then |
746 |
nInteractions = MnM_Map%interactionCount |
747 |
|
748 |
do i = 1, nInteractions |
749 |
MnM_Map%interactions(i)%shiftedPot = shiftedPot |
750 |
MnM_Map%interactions(i)%shiftedFrc = shiftedFrc |
751 |
MnM_Map%interactions(i)%rCut = thisRcut |
752 |
MnM_Map%interactions(i)%rCutWasSet = .true. |
753 |
enddo |
754 |
end if |
755 |
end if |
756 |
|
757 |
end subroutine setMnMDefaultCutoff |
758 |
|
759 |
subroutine copyAllData(v1, v2) |
760 |
type(MnMinteractionMap), pointer :: v1 |
761 |
type(MnMinteractionMap), pointer :: v2 |
762 |
integer :: i, j |
763 |
|
764 |
do i = 1, v1%interactionCount |
765 |
v2%interactions(i) = v1%interactions(i) |
766 |
enddo |
767 |
|
768 |
v2%interactionCount = v1%interactionCount |
769 |
return |
770 |
end subroutine copyAllData |
771 |
|
772 |
subroutine addInteraction(myInteraction) |
773 |
type(MNMtype) :: myInteraction |
774 |
type(MnMinteraction) :: nt |
775 |
integer :: id |
776 |
|
777 |
|
778 |
|
779 |
|
780 |
|
781 |
nt%interaction_type = myInteraction%MNMInteractionType |
782 |
nt%metal_atid = & |
783 |
getFirstMatchingElement(atypes, "c_ident", myInteraction%metal_atid) |
784 |
nt%nonmetal_atid = & |
785 |
getFirstMatchingElement(atypes, "c_ident", myInteraction%nonmetal_atid) |
786 |
|
787 |
|
788 |
select case (nt%interaction_type) |
789 |
case (MNM_LENNARDJONES) |
790 |
nt%sigma = myInteraction%sigma |
791 |
nt%epsilon = myInteraction%epsilon |
792 |
case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE) |
793 |
nt%R0 = myInteraction%R0 |
794 |
nt%D0 = myInteraction%D0 |
795 |
nt%beta0 = myInteraction%beta0 |
796 |
case(MNM_MAW) |
797 |
nt%R0 = myInteraction%R0 |
798 |
nt%D0 = myInteraction%D0 |
799 |
nt%beta0 = myInteraction%beta0 |
800 |
nt%alpha = myInteraction%alpha |
801 |
case default |
802 |
call handleError("MNM", "Unknown Interaction type") |
803 |
end select |
804 |
|
805 |
if (.not. associated(MnM_Map)) then |
806 |
call ensureCapacityHelper(MnM_Map, 1) |
807 |
else |
808 |
call ensureCapacityHelper(MnM_Map, MnM_Map%interactionCount + 1) |
809 |
end if |
810 |
|
811 |
MnM_Map%interactionCount = MnM_Map%interactionCount + 1 |
812 |
id = MnM_Map%interactionCount |
813 |
MnM_Map%interactions(id) = nt |
814 |
end subroutine addInteraction |
815 |
|
816 |
subroutine ensureCapacityHelper(this, minCapacity) |
817 |
type(MnMinteractionMap), pointer :: this, that |
818 |
integer, intent(in) :: minCapacity |
819 |
integer :: oldCapacity |
820 |
integer :: newCapacity |
821 |
logical :: resizeFlag |
822 |
|
823 |
resizeFlag = .false. |
824 |
|
825 |
! first time: allocate a new vector with default size |
826 |
|
827 |
if (.not. associated(this)) then |
828 |
this => MnMinitialize(minCapacity, 0) |
829 |
endif |
830 |
|
831 |
oldCapacity = size(this%interactions) |
832 |
|
833 |
if (minCapacity > oldCapacity) then |
834 |
if (this%capacityIncrement .gt. 0) then |
835 |
newCapacity = oldCapacity + this%capacityIncrement |
836 |
else |
837 |
newCapacity = oldCapacity * 2 |
838 |
endif |
839 |
if (newCapacity .lt. minCapacity) then |
840 |
newCapacity = minCapacity |
841 |
endif |
842 |
resizeFlag = .true. |
843 |
else |
844 |
newCapacity = oldCapacity |
845 |
endif |
846 |
|
847 |
if (resizeFlag) then |
848 |
that => MnMinitialize(newCapacity, this%capacityIncrement) |
849 |
call copyAllData(this, that) |
850 |
this => MnMdestroy(this) |
851 |
this => that |
852 |
endif |
853 |
end subroutine ensureCapacityHelper |
854 |
|
855 |
function MnMinitialize(cap, capinc) result(this) |
856 |
integer, intent(in) :: cap, capinc |
857 |
integer :: error |
858 |
type(MnMinteractionMap), pointer :: this |
859 |
|
860 |
nullify(this) |
861 |
|
862 |
if (cap < 0) then |
863 |
write(*,*) 'Bogus Capacity:', cap |
864 |
return |
865 |
endif |
866 |
allocate(this,stat=error) |
867 |
if ( error /= 0 ) then |
868 |
write(*,*) 'Could not allocate MnMinteractionMap!' |
869 |
return |
870 |
end if |
871 |
|
872 |
this%initialCapacity = cap |
873 |
this%capacityIncrement = capinc |
874 |
|
875 |
allocate(this%interactions(this%initialCapacity), stat=error) |
876 |
if(error /= 0) write(*,*) 'Could not allocate MnMinteraction!' |
877 |
|
878 |
end function MnMinitialize |
879 |
|
880 |
subroutine createInteractionLookup(this) |
881 |
type (MnMInteractionMap),pointer :: this |
882 |
integer :: biggestAtid, i, metal_atid, nonmetal_atid, error |
883 |
|
884 |
biggestAtid=-1 |
885 |
do i = 1, this%interactionCount |
886 |
metal_atid = this%interactions(i)%metal_atid |
887 |
nonmetal_atid = this%interactions(i)%nonmetal_atid |
888 |
|
889 |
if (metal_atid .gt. biggestAtid) biggestAtid = metal_atid |
890 |
if (nonmetal_atid .gt. biggestAtid) biggestAtid = nonmetal_atid |
891 |
enddo |
892 |
|
893 |
allocate(MnMinteractionLookup(biggestAtid,biggestAtid), stat=error) |
894 |
if (error /= 0) write(*,*) 'Could not allocate MnMinteractionLookup' |
895 |
|
896 |
do i = 1, this%interactionCount |
897 |
metal_atid = this%interactions(i)%metal_atid |
898 |
nonmetal_atid = this%interactions(i)%nonmetal_atid |
899 |
|
900 |
MnMinteractionLookup(metal_atid, nonmetal_atid) = i |
901 |
MnMinteractionLookup(nonmetal_atid, metal_atid) = i |
902 |
enddo |
903 |
end subroutine createInteractionLookup |
904 |
|
905 |
|
906 |
function MnMdestroy(this) result(null_this) |
907 |
logical :: done |
908 |
type(MnMinteractionMap), pointer :: this |
909 |
type(MnMinteractionMap), pointer :: null_this |
910 |
|
911 |
if (.not. associated(this)) then |
912 |
null_this => null() |
913 |
return |
914 |
end if |
915 |
|
916 |
!! Walk down the list and deallocate each of the map's components |
917 |
if(associated(this%interactions)) then |
918 |
deallocate(this%interactions) |
919 |
this%interactions=>null() |
920 |
endif |
921 |
deallocate(this) |
922 |
this => null() |
923 |
null_this => null() |
924 |
end function MnMdestroy |
925 |
|
926 |
|
927 |
subroutine deleteInteractions() |
928 |
MnM_Map => MnMdestroy(MnM_Map) |
929 |
return |
930 |
end subroutine deleteInteractions |
931 |
|
932 |
|
933 |
subroutine getLJfunc(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, ri12, ri13 |
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 |
ri12 = ri6*ri6 |
946 |
ri13 = ri12*ri |
947 |
|
948 |
myPot = 4.0_DP * (ri12 - ri6) |
949 |
myDeriv = 24.0_DP * (ri7 - 2.0_DP * ri13) |
950 |
|
951 |
return |
952 |
end subroutine getLJfunc |
953 |
|
954 |
subroutine getSoftFunc(r, myPot, myDeriv) |
955 |
|
956 |
real(kind=dp), intent(in) :: r |
957 |
real(kind=dp), intent(inout) :: myPot, myDeriv |
958 |
real(kind=dp) :: ri, ri2, ri6, ri7 |
959 |
real(kind=dp) :: a, b, c, d, dx |
960 |
integer :: j |
961 |
|
962 |
ri = 1.0_DP / r |
963 |
ri2 = ri*ri |
964 |
ri6 = ri2*ri2*ri2 |
965 |
ri7 = ri6*ri |
966 |
myPot = 4.0_DP * (ri6) |
967 |
myDeriv = - 24.0_DP * ri7 |
968 |
|
969 |
return |
970 |
end subroutine getSoftFunc |
971 |
|
972 |
|
973 |
|
974 |
|
975 |
|
976 |
|
977 |
end module MetalNonMetal |