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.2 2007-07-17 18:54:47 chuckv Exp $, $Date: 2007-07-17 18:54:47 $, $Name: not supported by cvs2svn $, $Revision: 1.2 $ |
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(inout) :: atom1, atom2 |
116 |
integer :: atid1, atid2, ljt1, ljt2 |
117 |
real( kind = dp ), intent(inout) :: rij, r2, rcut |
118 |
real( kind = dp ), intent(inout) :: pot, sw, vpair |
119 |
real( kind = dp ), intent(inout), dimension(3,nLocal) :: f |
120 |
real (kind=dp),intent(inout), dimension(9,nLocal) :: A |
121 |
real (kind=dp),intent(inout), dimension(3,nLocal) :: t |
122 |
real( kind = dp ), intent(inout), dimension(3) :: d |
123 |
real( kind = dp ), intent(inout), dimension(3) :: fpair |
124 |
logical, intent(inout) :: 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 |
interaction_id = MnMinteractionLookup(atid1, atid2) |
138 |
interaction_type = MnM_Map%interactions(interaction_id)%interaction_type |
139 |
|
140 |
select case (interaction_type) |
141 |
case (MNM_LENNARDJONES) |
142 |
call calc_mnm_lennardjones(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, & |
143 |
Pot, F, Do_pot, interaction_id) |
144 |
case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE) |
145 |
call calc_mnm_morse(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, & |
146 |
Pot, F, Do_pot, interaction_id,interaction_type) |
147 |
case(MNM_MAW) |
148 |
call calc_mnm_maw(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, & |
149 |
Pot,A, F,t, Do_pot, interaction_id) |
150 |
end select |
151 |
|
152 |
end subroutine do_mnm_pair |
153 |
|
154 |
subroutine calc_mnm_lennardjones(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, & |
155 |
Pot, F, Do_pot, interaction_id) |
156 |
|
157 |
integer, intent(inout) :: atom1, atom2 |
158 |
real( kind = dp ), intent(inout) :: rij, r2, rcut |
159 |
real( kind = dp ), intent(inout) :: pot, sw, vpair |
160 |
real( kind = dp ), intent(inout), dimension(3,nLocal) :: f |
161 |
real( kind = dp ), intent(inout), dimension(3) :: d |
162 |
real( kind = dp ), intent(inout), dimension(3) :: fpair |
163 |
logical, intent(inout) :: do_pot |
164 |
integer, intent(in) :: interaction_id |
165 |
|
166 |
! local Variables |
167 |
real( kind = dp ) :: drdx, drdy, drdz |
168 |
real( kind = dp ) :: fx, fy, fz |
169 |
real( kind = dp ) :: myPot, myPotC, myDeriv, myDerivC, ros, rcos |
170 |
real( kind = dp ) :: pot_temp, dudr |
171 |
real( kind = dp ) :: sigmai |
172 |
real( kind = dp ) :: epsilon |
173 |
logical :: isSoftCore, shiftedPot, shiftedFrc |
174 |
integer :: id1, id2, localError |
175 |
|
176 |
|
177 |
sigmai = MnM_Map%interactions(interaction_id)%sigma |
178 |
epsilon = MnM_Map%interactions(interaction_id)%epsilon |
179 |
shiftedPot = MnM_Map%interactions(interaction_id)%shiftedPot |
180 |
shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc |
181 |
|
182 |
ros = rij * sigmai |
183 |
|
184 |
|
185 |
call getLJfunc(ros, myPot, myDeriv) |
186 |
|
187 |
if (shiftedPot) then |
188 |
rcos = rcut * sigmai |
189 |
call getLJfunc(rcos, myPotC, myDerivC) |
190 |
myDerivC = 0.0_dp |
191 |
elseif (shiftedFrc) then |
192 |
rcos = rcut * sigmai |
193 |
call getLJfunc(rcos, myPotC, myDerivC) |
194 |
myPotC = myPotC + myDerivC * (rij - rcut) * sigmai |
195 |
else |
196 |
myPotC = 0.0_dp |
197 |
myDerivC = 0.0_dp |
198 |
endif |
199 |
|
200 |
|
201 |
|
202 |
pot_temp = epsilon * (myPot - myPotC) |
203 |
vpair = vpair + pot_temp |
204 |
dudr = sw * epsilon * (myDeriv - myDerivC) * sigmai |
205 |
|
206 |
drdx = d(1) / rij |
207 |
drdy = d(2) / rij |
208 |
drdz = d(3) / rij |
209 |
|
210 |
fx = dudr * drdx |
211 |
fy = dudr * drdy |
212 |
fz = dudr * drdz |
213 |
|
214 |
#ifdef IS_MPI |
215 |
if (do_pot) then |
216 |
pot_Row(VDW_POT,atom1) = pot_Row(VDW_POT,atom1) + sw*pot_temp*0.5 |
217 |
pot_Col(VDW_POT,atom2) = pot_Col(VDW_POT,atom2) + sw*pot_temp*0.5 |
218 |
endif |
219 |
|
220 |
f_Row(1,atom1) = f_Row(1,atom1) + fx |
221 |
f_Row(2,atom1) = f_Row(2,atom1) + fy |
222 |
f_Row(3,atom1) = f_Row(3,atom1) + fz |
223 |
|
224 |
f_Col(1,atom2) = f_Col(1,atom2) - fx |
225 |
f_Col(2,atom2) = f_Col(2,atom2) - fy |
226 |
f_Col(3,atom2) = f_Col(3,atom2) - fz |
227 |
|
228 |
#else |
229 |
if (do_pot) pot = pot + sw*pot_temp |
230 |
|
231 |
f(1,atom1) = f(1,atom1) + fx |
232 |
f(2,atom1) = f(2,atom1) + fy |
233 |
f(3,atom1) = f(3,atom1) + fz |
234 |
|
235 |
f(1,atom2) = f(1,atom2) - fx |
236 |
f(2,atom2) = f(2,atom2) - fy |
237 |
f(3,atom2) = f(3,atom2) - fz |
238 |
#endif |
239 |
|
240 |
#ifdef IS_MPI |
241 |
id1 = AtomRowToGlobal(atom1) |
242 |
id2 = AtomColToGlobal(atom2) |
243 |
#else |
244 |
id1 = atom1 |
245 |
id2 = atom2 |
246 |
#endif |
247 |
|
248 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
249 |
|
250 |
fpair(1) = fpair(1) + fx |
251 |
fpair(2) = fpair(2) + fy |
252 |
fpair(3) = fpair(3) + fz |
253 |
|
254 |
endif |
255 |
|
256 |
return |
257 |
|
258 |
|
259 |
end subroutine calc_mnm_lennardjones |
260 |
|
261 |
|
262 |
|
263 |
|
264 |
|
265 |
|
266 |
subroutine calc_mnm_morse(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, & |
267 |
Pot, f, Do_pot, interaction_id, interaction_type) |
268 |
integer, intent(inout) :: atom1, atom2 |
269 |
real( kind = dp ), intent(inout) :: rij, r2, rcut |
270 |
real( kind = dp ), intent(inout) :: pot, sw, vpair |
271 |
real( kind = dp ), intent(inout), dimension(3,nLocal) :: f |
272 |
real( kind = dp ), intent(inout), dimension(3) :: d |
273 |
real( kind = dp ), intent(inout), dimension(3) :: fpair |
274 |
logical, intent(inout) :: do_pot |
275 |
integer, intent(in) :: interaction_id, interaction_type |
276 |
|
277 |
! Local Variables |
278 |
real(kind=dp) :: Beta0 |
279 |
real(kind=dp) :: R0 |
280 |
real(kind=dp) :: D0 |
281 |
real(kind=dp) :: expt |
282 |
real(kind=dp) :: expt2 |
283 |
real(kind=dp) :: expfnc |
284 |
real(kind=dp) :: expfnc2 |
285 |
real(kind=dp) :: D_expt |
286 |
real(kind=dp) :: D_expt2 |
287 |
real(kind=dp) :: rcos |
288 |
real(kind=dp) :: exptC |
289 |
real(kind=dp) :: expt2C |
290 |
real(kind=dp) :: expfncC |
291 |
real(kind=dp) :: expfnc2C |
292 |
real(kind=dp) :: D_expt2C |
293 |
real(kind=dp) :: D_exptC |
294 |
|
295 |
real(kind=dp) :: myPot |
296 |
real(kind=dp) :: myPotC |
297 |
real(kind=dp) :: myDeriv |
298 |
real(kind=dp) :: myDerivC |
299 |
real(kind=dp) :: pot_temp |
300 |
real(kind=dp) :: fx,fy,fz |
301 |
real(kind=dp) :: drdx,drdy,drdz |
302 |
real(kind=dp) :: dudr |
303 |
integer :: id1,id2 |
304 |
|
305 |
|
306 |
D0 = MnM_Map%interactions(interaction_id)%D0 |
307 |
R0 = MnM_Map%interactions(interaction_id)%r0 |
308 |
Beta0 = MnM_Map%interactions(interaction_id)%Beta0 |
309 |
! shiftedFrc = MnM_Map%interactions(interaction_id)%shiftedFrc |
310 |
|
311 |
|
312 |
! V(r) = D_e exp(-a(r-re)(exp(-a(r-re))-2) |
313 |
|
314 |
expt = -R0*(rij-R0) |
315 |
expt2 = 2.0_dp*expt |
316 |
expfnc = exp(expt) |
317 |
expfnc2 = exp(expt2) |
318 |
D_expt = D0*expt |
319 |
D_expt2 = D0*expt2 |
320 |
|
321 |
rcos = rcut*Beta0 |
322 |
exptC = -R0*(rcos-R0) |
323 |
expt2C = 2.0_dp*expt |
324 |
expfncC = exp(expt) |
325 |
expfnc2C = exp(expt2) |
326 |
D_exptC = D0*expt |
327 |
D_expt2C = D0*expt2 |
328 |
|
329 |
select case (interaction_type) |
330 |
|
331 |
|
332 |
case (MNM_SHIFTEDMORSE) |
333 |
|
334 |
myPot = D_expt * (D_expt - 2.0_dp) |
335 |
myPotC = D_exptC * (D_exptC - 2.0_dp) |
336 |
|
337 |
myDeriv = -(D_expt2 - D_expt) * (-2.0_dp + D_expt) |
338 |
myDerivC = -(D_expt2C - D_exptC) * (-2.0_dp + D_exptC) |
339 |
|
340 |
case (MNM_REPULSIVEMORSE) |
341 |
|
342 |
myPot = D_expt2 |
343 |
myPotC = D_expt2C |
344 |
|
345 |
myDeriv = -2.0_dp * D_expt2 |
346 |
myDerivC = -2.0_dp * D_expt2C |
347 |
|
348 |
end select |
349 |
|
350 |
myPotC = myPotC + myDerivC*(rij - rcut)*Beta0 |
351 |
|
352 |
pot_temp = (myPot - myPotC) |
353 |
vpair = vpair + pot_temp |
354 |
dudr = sw * (myDeriv - myDerivC) |
355 |
|
356 |
drdx = d(1) / rij |
357 |
drdy = d(2) / rij |
358 |
drdz = d(3) / rij |
359 |
|
360 |
fx = dudr * drdx |
361 |
fy = dudr * drdy |
362 |
fz = dudr * drdz |
363 |
|
364 |
#ifdef IS_MPI |
365 |
if (do_pot) then |
366 |
pot_Row(VDW_POT,atom1) = pot_Row(VDW_POT,atom1) + sw*pot_temp*0.5 |
367 |
pot_Col(VDW_POT,atom2) = pot_Col(VDW_POT,atom2) + sw*pot_temp*0.5 |
368 |
endif |
369 |
|
370 |
f_Row(1,atom1) = f_Row(1,atom1) + fx |
371 |
f_Row(2,atom1) = f_Row(2,atom1) + fy |
372 |
f_Row(3,atom1) = f_Row(3,atom1) + fz |
373 |
|
374 |
f_Col(1,atom2) = f_Col(1,atom2) - fx |
375 |
f_Col(2,atom2) = f_Col(2,atom2) - fy |
376 |
f_Col(3,atom2) = f_Col(3,atom2) - fz |
377 |
|
378 |
#else |
379 |
if (do_pot) pot = pot + sw*pot_temp |
380 |
|
381 |
f(1,atom1) = f(1,atom1) + fx |
382 |
f(2,atom1) = f(2,atom1) + fy |
383 |
f(3,atom1) = f(3,atom1) + fz |
384 |
|
385 |
f(1,atom2) = f(1,atom2) - fx |
386 |
f(2,atom2) = f(2,atom2) - fy |
387 |
f(3,atom2) = f(3,atom2) - fz |
388 |
#endif |
389 |
|
390 |
#ifdef IS_MPI |
391 |
id1 = AtomRowToGlobal(atom1) |
392 |
id2 = AtomColToGlobal(atom2) |
393 |
#else |
394 |
id1 = atom1 |
395 |
id2 = atom2 |
396 |
#endif |
397 |
|
398 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
399 |
|
400 |
fpair(1) = fpair(1) + fx |
401 |
fpair(2) = fpair(2) + fy |
402 |
fpair(3) = fpair(3) + fz |
403 |
|
404 |
endif |
405 |
|
406 |
return |
407 |
|
408 |
end subroutine calc_mnm_morse |
409 |
|
410 |
|
411 |
|
412 |
|
413 |
subroutine calc_mnm_maw(Atom1, Atom2, D, Rij, R2, Rcut, Sw, Vpair, Fpair, & |
414 |
Pot,A, F,t, Do_pot, interaction_id) |
415 |
integer, intent(inout) :: atom1, atom2 |
416 |
real( kind = dp ), intent(inout) :: rij, r2, rcut |
417 |
real( kind = dp ), intent(inout) :: pot, sw, vpair |
418 |
real( kind = dp ), intent(inout), dimension(3,nLocal) :: f |
419 |
real (kind=dp),intent(inout), dimension(9,nLocal) :: A |
420 |
real (kind=dp),intent(inout), dimension(3,nLocal) :: t |
421 |
|
422 |
real( kind = dp ), intent(inout), dimension(3) :: d |
423 |
real( kind = dp ), intent(inout), dimension(3) :: fpair |
424 |
logical, intent(inout) :: do_pot |
425 |
|
426 |
integer, intent(in) :: interaction_id |
427 |
|
428 |
end subroutine calc_mnm_maw |
429 |
|
430 |
|
431 |
subroutine setMnMDefaultCutoff(thisRcut, shiftedPot, shiftedFrc) |
432 |
real(kind=dp), intent(in) :: thisRcut |
433 |
logical, intent(in) :: shiftedPot |
434 |
logical, intent(in) :: shiftedFrc |
435 |
integer i, nInteractions |
436 |
defaultCutoff = thisRcut |
437 |
defaultShiftPot = shiftedPot |
438 |
defaultShiftFrc = shiftedFrc |
439 |
|
440 |
if(MnM_Map%interactionCount /= 0) then |
441 |
nInteractions = MnM_Map%interactionCount |
442 |
|
443 |
do i = 1, nInteractions |
444 |
MnM_Map%interactions(i)%shiftedPot = shiftedPot |
445 |
MnM_Map%interactions(i)%shiftedFrc = shiftedFrc |
446 |
MnM_Map%interactions(i)%rCut = thisRcut |
447 |
MnM_Map%interactions(i)%rCutWasSet = .true. |
448 |
enddo |
449 |
end if |
450 |
|
451 |
end subroutine setMnMDefaultCutoff |
452 |
|
453 |
subroutine copyAllData(v1, v2) |
454 |
type(MnMinteractionMap), pointer :: v1 |
455 |
type(MnMinteractionMap), pointer :: v2 |
456 |
integer :: i, j |
457 |
|
458 |
do i = 1, v1%interactionCount |
459 |
v2%interactions(i) = v1%interactions(i) |
460 |
enddo |
461 |
|
462 |
v2%interactionCount = v1%interactionCount |
463 |
return |
464 |
end subroutine copyAllData |
465 |
|
466 |
subroutine addInteraction(myInteraction) |
467 |
type(MNMtype) :: myInteraction |
468 |
type(MnMinteraction) :: nt |
469 |
integer :: id |
470 |
|
471 |
nt%interaction_type = myInteraction%MNMInteractionType |
472 |
nt%metal_atid = myInteraction%metal_atid |
473 |
nt%nonmetal_atid = myInteraction%nonmetal_atid |
474 |
|
475 |
select case (nt%interaction_type) |
476 |
case (MNM_LENNARDJONES) |
477 |
nt%sigma = myInteraction%sigma |
478 |
nt%epsilon = myInteraction%epsilon |
479 |
case(MNM_REPULSIVEMORSE, MNM_SHIFTEDMORSE) |
480 |
nt%R0 = myInteraction%R0 |
481 |
nt%D0 = myInteraction%D0 |
482 |
nt%beta0 = myInteraction%beta0 |
483 |
case(MNM_MAW) |
484 |
nt%R0 = myInteraction%R0 |
485 |
nt%D0 = myInteraction%D0 |
486 |
nt%beta0 = myInteraction%beta0 |
487 |
nt%betaH = myInteraction%betaH |
488 |
nt%alpha = myInteraction%alpha |
489 |
nt%gamma = myInteraction%gamma |
490 |
case default |
491 |
call handleError("MNM", "Unknown Interaction type") |
492 |
end select |
493 |
|
494 |
if (.not. associated(MnM_Map)) then |
495 |
call ensureCapacityHelper(MnM_Map, 1) |
496 |
else |
497 |
call ensureCapacityHelper(MnM_Map, MnM_Map%interactionCount + 1) |
498 |
end if |
499 |
|
500 |
MnM_Map%interactionCount = MnM_Map%interactionCount + 1 |
501 |
id = MnM_Map%interactionCount |
502 |
MnM_Map%interactions(id) = nt |
503 |
end subroutine addInteraction |
504 |
|
505 |
subroutine ensureCapacityHelper(this, minCapacity) |
506 |
type(MnMinteractionMap), pointer :: this, that |
507 |
integer, intent(in) :: minCapacity |
508 |
integer :: oldCapacity |
509 |
integer :: newCapacity |
510 |
logical :: resizeFlag |
511 |
|
512 |
resizeFlag = .false. |
513 |
|
514 |
! first time: allocate a new vector with default size |
515 |
|
516 |
if (.not. associated(this)) then |
517 |
this => MnMinitialize(minCapacity, 0) |
518 |
endif |
519 |
|
520 |
oldCapacity = size(this%interactions) |
521 |
|
522 |
if (minCapacity > oldCapacity) then |
523 |
if (this%capacityIncrement .gt. 0) then |
524 |
newCapacity = oldCapacity + this%capacityIncrement |
525 |
else |
526 |
newCapacity = oldCapacity * 2 |
527 |
endif |
528 |
if (newCapacity .lt. minCapacity) then |
529 |
newCapacity = minCapacity |
530 |
endif |
531 |
resizeFlag = .true. |
532 |
else |
533 |
newCapacity = oldCapacity |
534 |
endif |
535 |
|
536 |
if (resizeFlag) then |
537 |
that => MnMinitialize(newCapacity, this%capacityIncrement) |
538 |
call copyAllData(this, that) |
539 |
this => MnMdestroy(this) |
540 |
this => that |
541 |
endif |
542 |
end subroutine ensureCapacityHelper |
543 |
|
544 |
function MnMinitialize(cap, capinc) result(this) |
545 |
integer, intent(in) :: cap, capinc |
546 |
integer :: error |
547 |
type(MnMinteractionMap), pointer :: this |
548 |
|
549 |
nullify(this) |
550 |
|
551 |
if (cap < 0) then |
552 |
write(*,*) 'Bogus Capacity:', cap |
553 |
return |
554 |
endif |
555 |
allocate(this,stat=error) |
556 |
if ( error /= 0 ) then |
557 |
write(*,*) 'Could not allocate MnMinteractionMap!' |
558 |
return |
559 |
end if |
560 |
|
561 |
this%initialCapacity = cap |
562 |
this%capacityIncrement = capinc |
563 |
|
564 |
allocate(this%interactions(this%initialCapacity), stat=error) |
565 |
if(error /= 0) write(*,*) 'Could not allocate MnMinteraction!' |
566 |
|
567 |
end function MnMinitialize |
568 |
|
569 |
subroutine createInteractionLookup(this) |
570 |
type(MnMinteractionMap), pointer :: this |
571 |
integer :: biggestAtid, i, metal_atid, nonmetal_atid, error |
572 |
|
573 |
biggestAtid=-1 |
574 |
do i = 1, this%interactionCount |
575 |
metal_atid = this%interactions(i)%metal_atid |
576 |
nonmetal_atid = this%interactions(i)%nonmetal_atid |
577 |
|
578 |
if (metal_atid .gt. biggestAtid) biggestAtid = metal_atid |
579 |
if (nonmetal_atid .gt. biggestAtid) biggestAtid = nonmetal_atid |
580 |
enddo |
581 |
|
582 |
allocate(MnMinteractionLookup(biggestAtid,biggestAtid), stat=error) |
583 |
if (error /= 0) write(*,*) 'Could not allocate MnMinteractionLookup' |
584 |
|
585 |
do i = 1, this%interactionCount |
586 |
metal_atid = this%interactions(i)%metal_atid |
587 |
nonmetal_atid = this%interactions(i)%nonmetal_atid |
588 |
|
589 |
MnMinteractionLookup(metal_atid, nonmetal_atid) = i |
590 |
MnMinteractionLookup(nonmetal_atid, metal_atid) = i |
591 |
enddo |
592 |
end subroutine createInteractionLookup |
593 |
|
594 |
|
595 |
function MnMdestroy(this) result(null_this) |
596 |
logical :: done |
597 |
type(MnMinteractionMap), pointer :: this |
598 |
type(MnMinteractionMap), pointer :: null_this |
599 |
|
600 |
if (.not. associated(this)) then |
601 |
null_this => null() |
602 |
return |
603 |
end if |
604 |
|
605 |
!! Walk down the list and deallocate each of the map's components |
606 |
if(associated(this%interactions)) then |
607 |
deallocate(this%interactions) |
608 |
this%interactions=>null() |
609 |
endif |
610 |
deallocate(this) |
611 |
this => null() |
612 |
null_this => null() |
613 |
end function MnMdestroy |
614 |
|
615 |
|
616 |
subroutine deleteInteractions() |
617 |
MnM_Map => MnMdestroy(MnM_Map) |
618 |
return |
619 |
end subroutine deleteInteractions |
620 |
|
621 |
|
622 |
subroutine getLJfunc(r, myPot, myDeriv) |
623 |
|
624 |
real(kind=dp), intent(in) :: r |
625 |
real(kind=dp), intent(inout) :: myPot, myDeriv |
626 |
real(kind=dp) :: ri, ri2, ri6, ri7, ri12, ri13 |
627 |
real(kind=dp) :: a, b, c, d, dx |
628 |
integer :: j |
629 |
|
630 |
ri = 1.0_DP / r |
631 |
ri2 = ri*ri |
632 |
ri6 = ri2*ri2*ri2 |
633 |
ri7 = ri6*ri |
634 |
ri12 = ri6*ri6 |
635 |
ri13 = ri12*ri |
636 |
|
637 |
myPot = 4.0_DP * (ri12 - ri6) |
638 |
myDeriv = 24.0_DP * (ri7 - 2.0_DP * ri13) |
639 |
|
640 |
return |
641 |
end subroutine getLJfunc |
642 |
|
643 |
subroutine getSoftFunc(r, myPot, myDeriv) |
644 |
|
645 |
real(kind=dp), intent(in) :: r |
646 |
real(kind=dp), intent(inout) :: myPot, myDeriv |
647 |
real(kind=dp) :: ri, ri2, ri6, ri7 |
648 |
real(kind=dp) :: a, b, c, d, dx |
649 |
integer :: j |
650 |
|
651 |
ri = 1.0_DP / r |
652 |
ri2 = ri*ri |
653 |
ri6 = ri2*ri2*ri2 |
654 |
ri7 = ri6*ri |
655 |
myPot = 4.0_DP * (ri6) |
656 |
myDeriv = - 24.0_DP * ri7 |
657 |
|
658 |
return |
659 |
end subroutine getSoftFunc |
660 |
|
661 |
|
662 |
|
663 |
|
664 |
|
665 |
|
666 |
end module MetalNonMetal |