1 |
!! |
2 |
!! Copyright (c) 2005, 2009 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. Redistributions of source code must retain the above copyright |
10 |
!! notice, this list of conditions and the following disclaimer. |
11 |
!! |
12 |
!! 2. Redistributions in binary form must reproduce the above copyright |
13 |
!! notice, this list of conditions and the following disclaimer in the |
14 |
!! documentation and/or other materials provided with the |
15 |
!! distribution. |
16 |
!! |
17 |
!! This software is provided "AS IS," without a warranty of any |
18 |
!! kind. All express or implied conditions, representations and |
19 |
!! warranties, including any implied warranty of merchantability, |
20 |
!! fitness for a particular purpose or non-infringement, are hereby |
21 |
!! excluded. The University of Notre Dame and its licensors shall not |
22 |
!! be liable for any damages suffered by licensee as a result of |
23 |
!! using, modifying or distributing the software or its |
24 |
!! derivatives. In no event will the University of Notre Dame or its |
25 |
!! licensors be liable for any lost revenue, profit or data, or for |
26 |
!! direct, indirect, special, consequential, incidental or punitive |
27 |
!! damages, however caused and regardless of the theory of liability, |
28 |
!! arising out of the use of or inability to use software, even if the |
29 |
!! University of Notre Dame has been advised of the possibility of |
30 |
!! such damages. |
31 |
!! |
32 |
!! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
!! research, please cite the appropriate papers when you publish your |
34 |
!! work. Good starting points are: |
35 |
!! |
36 |
!! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
!! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
!! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). |
39 |
!! [4] Vardeman & Gezelter, in progress (2009). |
40 |
!! |
41 |
|
42 |
!! doForces.F90 |
43 |
!! module doForces |
44 |
!! Calculates Long Range forces. |
45 |
|
46 |
!! @author Charles F. Vardeman II |
47 |
!! @author Matthew Meineke |
48 |
!! @version $Id$, $Date$, $Name: not supported by cvs2svn $, $Revision$ |
49 |
|
50 |
|
51 |
module doForces |
52 |
use force_globals |
53 |
use fForceOptions |
54 |
use simulation |
55 |
use definitions |
56 |
use atype_module |
57 |
use switcheroo |
58 |
use neighborLists |
59 |
use sticky |
60 |
use electrostatic_module |
61 |
use gayberne |
62 |
use shapes |
63 |
use vector_class |
64 |
use eam |
65 |
use MetalNonMetal |
66 |
use suttonchen |
67 |
use status |
68 |
use ISO_C_BINDING |
69 |
|
70 |
#ifdef IS_MPI |
71 |
use mpiSimulation |
72 |
#endif |
73 |
|
74 |
implicit none |
75 |
PRIVATE |
76 |
|
77 |
!!$ interface |
78 |
!!$ function getSigma(atid) |
79 |
!!$ import :: c_int |
80 |
!!$ import :: c_double |
81 |
!!$ integer(c_int), intent(in) :: atid |
82 |
!!$ real(c_double) :: getSigma |
83 |
!!$ end function getSigma |
84 |
!!$ function getEpsilon(atid) |
85 |
!!$ import :: c_int |
86 |
!!$ import :: c_double |
87 |
!!$ integer(c_int), intent(in) :: atid |
88 |
!!$ real(c_double) :: getEpsilon |
89 |
!!$ end function getEpsilon |
90 |
!!$ subroutine do_lj_pair(atid1, atid2, d, rij, r2, rcut, sw, vdwMult, & |
91 |
!!$ vpair, pot, f1) |
92 |
!!$ import :: c_int |
93 |
!!$ import :: c_double |
94 |
!!$ integer(c_int), intent(in) :: atid1, atid2 |
95 |
!!$ real(c_double), intent(in) :: rij, r2, rcut, vdwMult, sw |
96 |
!!$ real(c_double), intent(in), dimension(3) :: d |
97 |
!!$ real(c_double), intent(inout) :: pot, vpair |
98 |
!!$ real(c_double), intent(inout), dimension(3) :: f1 |
99 |
!!$ end subroutine do_lj_pair |
100 |
!!$ end interface |
101 |
!!$ interface |
102 |
!!$ function getSigmaFunc() bind(c,name="getSigma") |
103 |
!!$ import :: c_funptr |
104 |
!!$ type(c_funptr) :: getSigmaFunc |
105 |
!!$ end function getSigmaFunc |
106 |
!!$ function getEpsilonFunc() bind(c,name="getEpsilon") |
107 |
!!$ import :: c_funptr |
108 |
!!$ type(c_funptr) :: getEpsilonFunc |
109 |
!!$ end function getEpsilonFunc |
110 |
!!$ function getLJPfunc() bind(c,name="LJ_do_pair") |
111 |
!!$ import :: c_funptr |
112 |
!!$ type(c_funptr) :: getLJPfunc |
113 |
!!$ end function getLJPfunc |
114 |
!!$ end interface |
115 |
!!$ |
116 |
!!$ type (c_funptr) :: gsfunptr, gefunptr, dljpsubptr |
117 |
!!$ procedure(func), pointer :: getSigma, getEpsilon, do_lj_pair |
118 |
|
119 |
real(kind=dp), external :: getSigma |
120 |
real(kind=dp), external :: getEpsilon |
121 |
|
122 |
#define __FORTRAN90 |
123 |
#include "UseTheForce/fCutoffPolicy.h" |
124 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
125 |
#include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h" |
126 |
|
127 |
INTEGER, PARAMETER:: PREPAIR_LOOP = 1 |
128 |
INTEGER, PARAMETER:: PAIR_LOOP = 2 |
129 |
|
130 |
logical, save :: haveNeighborList = .false. |
131 |
logical, save :: haveSIMvariables = .false. |
132 |
logical, save :: haveSaneForceField = .false. |
133 |
logical, save :: haveInteractionHash = .false. |
134 |
logical, save :: haveGtypeCutoffMap = .false. |
135 |
logical, save :: haveDefaultCutoffs = .false. |
136 |
logical, save :: haveSkinThickness = .false. |
137 |
logical, save :: haveElectrostaticSummationMethod = .false. |
138 |
logical, save :: haveCutoffPolicy = .false. |
139 |
logical, save :: VisitCutoffsAfterComputing = .false. |
140 |
logical, save :: do_box_dipole = .false. |
141 |
|
142 |
logical, save :: FF_uses_DirectionalAtoms |
143 |
logical, save :: FF_uses_Dipoles |
144 |
logical, save :: FF_uses_GayBerne |
145 |
logical, save :: FF_uses_EAM |
146 |
logical, save :: FF_uses_SC |
147 |
logical, save :: FF_uses_MNM |
148 |
|
149 |
|
150 |
logical, save :: SIM_uses_DirectionalAtoms |
151 |
logical, save :: SIM_uses_EAM |
152 |
logical, save :: SIM_uses_SC |
153 |
logical, save :: SIM_uses_MNM |
154 |
logical, save :: SIM_requires_postpair_calc |
155 |
logical, save :: SIM_requires_prepair_calc |
156 |
logical, save :: SIM_uses_PBC |
157 |
logical, save :: SIM_uses_AtomicVirial |
158 |
|
159 |
integer, save :: electrostaticSummationMethod |
160 |
integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY |
161 |
|
162 |
real(kind=dp), save :: defaultRcut, defaultRsw, largestRcut |
163 |
real(kind=dp), save :: skinThickness |
164 |
logical, save :: defaultDoShiftPot |
165 |
logical, save :: defaultDoShiftFrc |
166 |
|
167 |
public :: init_FF |
168 |
public :: setCutoffs |
169 |
public :: cWasLame |
170 |
public :: setElectrostaticMethod |
171 |
public :: setBoxDipole |
172 |
public :: getBoxDipole |
173 |
public :: setCutoffPolicy |
174 |
public :: setSkinThickness |
175 |
public :: do_force_loop |
176 |
|
177 |
#ifdef PROFILE |
178 |
public :: getforcetime |
179 |
real, save :: forceTime = 0 |
180 |
real :: forceTimeInitial, forceTimeFinal |
181 |
integer :: nLoops |
182 |
#endif |
183 |
|
184 |
!! Variables for cutoff mapping and interaction mapping |
185 |
! Bit hash to determine pair-pair interactions. |
186 |
integer, dimension(:,:), allocatable :: InteractionHash |
187 |
real(kind=dp), dimension(:), allocatable :: atypeMaxCutoff |
188 |
real(kind=dp), dimension(:), allocatable, target :: groupMaxCutoffRow |
189 |
real(kind=dp), dimension(:), pointer :: groupMaxCutoffCol |
190 |
|
191 |
integer, dimension(:), allocatable, target :: groupToGtypeRow |
192 |
integer, dimension(:), pointer :: groupToGtypeCol => null() |
193 |
|
194 |
real(kind=dp), dimension(:), allocatable,target :: gtypeMaxCutoffRow |
195 |
real(kind=dp), dimension(:), pointer :: gtypeMaxCutoffCol |
196 |
type ::gtypeCutoffs |
197 |
real(kind=dp) :: rcut |
198 |
real(kind=dp) :: rcutsq |
199 |
real(kind=dp) :: rlistsq |
200 |
end type gtypeCutoffs |
201 |
type(gtypeCutoffs), dimension(:,:), allocatable :: gtypeCutoffMap |
202 |
|
203 |
real(kind=dp), dimension(3) :: boxDipole |
204 |
|
205 |
contains |
206 |
|
207 |
subroutine createInteractionHash() |
208 |
integer :: nAtypes |
209 |
integer :: i |
210 |
integer :: j |
211 |
integer :: iHash |
212 |
!! Test Types |
213 |
logical :: i_is_LJ |
214 |
logical :: i_is_Elect |
215 |
logical :: i_is_Sticky |
216 |
logical :: i_is_StickyP |
217 |
logical :: i_is_GB |
218 |
logical :: i_is_EAM |
219 |
logical :: i_is_Shape |
220 |
logical :: i_is_SC |
221 |
logical :: j_is_LJ |
222 |
logical :: j_is_Elect |
223 |
logical :: j_is_Sticky |
224 |
logical :: j_is_StickyP |
225 |
logical :: j_is_GB |
226 |
logical :: j_is_EAM |
227 |
logical :: j_is_Shape |
228 |
logical :: j_is_SC |
229 |
real(kind=dp) :: myRcut |
230 |
|
231 |
if (.not. associated(atypes)) then |
232 |
call handleError("doForces", "atypes was not present before call of createInteractionHash!") |
233 |
return |
234 |
endif |
235 |
|
236 |
nAtypes = getSize(atypes) |
237 |
|
238 |
if (nAtypes == 0) then |
239 |
call handleError("doForces", "nAtypes was zero during call of createInteractionHash!") |
240 |
return |
241 |
end if |
242 |
|
243 |
if (.not. allocated(InteractionHash)) then |
244 |
allocate(InteractionHash(nAtypes,nAtypes)) |
245 |
else |
246 |
deallocate(InteractionHash) |
247 |
allocate(InteractionHash(nAtypes,nAtypes)) |
248 |
endif |
249 |
|
250 |
if (.not. allocated(atypeMaxCutoff)) then |
251 |
allocate(atypeMaxCutoff(nAtypes)) |
252 |
else |
253 |
deallocate(atypeMaxCutoff) |
254 |
allocate(atypeMaxCutoff(nAtypes)) |
255 |
endif |
256 |
|
257 |
do i = 1, nAtypes |
258 |
call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ) |
259 |
call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect) |
260 |
call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky) |
261 |
call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP) |
262 |
call getElementProperty(atypes, i, "is_GayBerne", i_is_GB) |
263 |
call getElementProperty(atypes, i, "is_EAM", i_is_EAM) |
264 |
call getElementProperty(atypes, i, "is_Shape", i_is_Shape) |
265 |
call getElementProperty(atypes, i, "is_SC", i_is_SC) |
266 |
|
267 |
do j = i, nAtypes |
268 |
|
269 |
iHash = 0 |
270 |
myRcut = 0.0_dp |
271 |
|
272 |
call getElementProperty(atypes, j, "is_LennardJones", j_is_LJ) |
273 |
call getElementProperty(atypes, j, "is_Electrostatic", j_is_Elect) |
274 |
call getElementProperty(atypes, j, "is_Sticky", j_is_Sticky) |
275 |
call getElementProperty(atypes, j, "is_StickyPower", j_is_StickyP) |
276 |
call getElementProperty(atypes, j, "is_GayBerne", j_is_GB) |
277 |
call getElementProperty(atypes, j, "is_EAM", j_is_EAM) |
278 |
call getElementProperty(atypes, j, "is_Shape", j_is_Shape) |
279 |
call getElementProperty(atypes, j, "is_SC", j_is_SC) |
280 |
|
281 |
if (i_is_LJ .and. j_is_LJ) then |
282 |
iHash = ior(iHash, LJ_PAIR) |
283 |
endif |
284 |
|
285 |
if (i_is_Elect .and. j_is_Elect) then |
286 |
iHash = ior(iHash, ELECTROSTATIC_PAIR) |
287 |
endif |
288 |
|
289 |
if (i_is_Sticky .and. j_is_Sticky) then |
290 |
iHash = ior(iHash, STICKY_PAIR) |
291 |
endif |
292 |
|
293 |
if (i_is_StickyP .and. j_is_StickyP) then |
294 |
iHash = ior(iHash, STICKYPOWER_PAIR) |
295 |
endif |
296 |
|
297 |
if (i_is_EAM .and. j_is_EAM) then |
298 |
iHash = ior(iHash, EAM_PAIR) |
299 |
endif |
300 |
|
301 |
if (i_is_SC .and. j_is_SC) then |
302 |
iHash = ior(iHash, SC_PAIR) |
303 |
endif |
304 |
|
305 |
if (i_is_GB .and. j_is_GB) iHash = ior(iHash, GAYBERNE_PAIR) |
306 |
if (i_is_GB .and. j_is_LJ) iHash = ior(iHash, GAYBERNE_LJ) |
307 |
if (i_is_LJ .and. j_is_GB) iHash = ior(iHash, GAYBERNE_LJ) |
308 |
|
309 |
if ((i_is_EAM.or.i_is_SC).and.(.not.(j_is_EAM.or.j_is_SC))) iHash = ior(iHash, MNM_PAIR) |
310 |
if ((j_is_EAM.or.j_is_SC).and.(.not.(i_is_EAM.or.i_is_SC))) iHash = ior(iHash, MNM_PAIR) |
311 |
|
312 |
if (i_is_Shape .and. j_is_Shape) iHash = ior(iHash, SHAPE_PAIR) |
313 |
if (i_is_Shape .and. j_is_LJ) iHash = ior(iHash, SHAPE_LJ) |
314 |
if (i_is_LJ .and. j_is_Shape) iHash = ior(iHash, SHAPE_LJ) |
315 |
|
316 |
|
317 |
InteractionHash(i,j) = iHash |
318 |
InteractionHash(j,i) = iHash |
319 |
|
320 |
end do |
321 |
|
322 |
end do |
323 |
|
324 |
haveInteractionHash = .true. |
325 |
end subroutine createInteractionHash |
326 |
|
327 |
subroutine createGtypeCutoffMap() |
328 |
|
329 |
logical :: i_is_LJ |
330 |
logical :: i_is_Elect |
331 |
logical :: i_is_Sticky |
332 |
logical :: i_is_StickyP |
333 |
logical :: i_is_GB |
334 |
logical :: i_is_EAM |
335 |
logical :: i_is_Shape |
336 |
logical :: i_is_SC |
337 |
logical :: GtypeFound |
338 |
|
339 |
integer :: myStatus, nAtypes, i, j, istart, iend, jstart, jend |
340 |
integer :: n_in_i, me_i, ia, g, atom1, ja, n_in_j,me_j |
341 |
integer :: nGroupsInRow |
342 |
integer :: nGroupsInCol |
343 |
integer :: nGroupTypesRow,nGroupTypesCol |
344 |
real(kind=dp):: thisSigma, bigSigma, thisRcut, tradRcut, tol |
345 |
real(kind=dp) :: biggestAtypeCutoff |
346 |
|
347 |
if (.not. haveInteractionHash) then |
348 |
call createInteractionHash() |
349 |
endif |
350 |
#ifdef IS_MPI |
351 |
nGroupsInRow = getNgroupsInRow(plan_group_row) |
352 |
nGroupsInCol = getNgroupsInCol(plan_group_col) |
353 |
#endif |
354 |
nAtypes = getSize(atypes) |
355 |
! Set all of the initial cutoffs to zero. |
356 |
atypeMaxCutoff = 0.0_dp |
357 |
biggestAtypeCutoff = 0.0_dp |
358 |
do i = 1, nAtypes |
359 |
if (SimHasAtype(i)) then |
360 |
call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ) |
361 |
call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect) |
362 |
call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky) |
363 |
call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP) |
364 |
call getElementProperty(atypes, i, "is_GayBerne", i_is_GB) |
365 |
call getElementProperty(atypes, i, "is_EAM", i_is_EAM) |
366 |
call getElementProperty(atypes, i, "is_Shape", i_is_Shape) |
367 |
call getElementProperty(atypes, i, "is_SC", i_is_SC) |
368 |
|
369 |
if (haveDefaultCutoffs) then |
370 |
atypeMaxCutoff(i) = defaultRcut |
371 |
else |
372 |
if (i_is_LJ) then |
373 |
thisRcut = getSigma(i) * 2.5_dp |
374 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
375 |
endif |
376 |
if (i_is_Elect) then |
377 |
thisRcut = defaultRcut |
378 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
379 |
endif |
380 |
if (i_is_Sticky) then |
381 |
thisRcut = getStickyCut(i) |
382 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
383 |
endif |
384 |
if (i_is_StickyP) then |
385 |
thisRcut = getStickyPowerCut(i) |
386 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
387 |
endif |
388 |
if (i_is_GB) then |
389 |
thisRcut = getGayBerneCut(i) |
390 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
391 |
endif |
392 |
if (i_is_EAM) then |
393 |
thisRcut = getEAMCut(i) |
394 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
395 |
endif |
396 |
if (i_is_Shape) then |
397 |
thisRcut = getShapeCut(i) |
398 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
399 |
endif |
400 |
if (i_is_SC) then |
401 |
thisRcut = getSCCut(i) |
402 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
403 |
endif |
404 |
endif |
405 |
|
406 |
if (atypeMaxCutoff(i).gt.biggestAtypeCutoff) then |
407 |
biggestAtypeCutoff = atypeMaxCutoff(i) |
408 |
endif |
409 |
|
410 |
endif |
411 |
enddo |
412 |
|
413 |
istart = 1 |
414 |
jstart = 1 |
415 |
#ifdef IS_MPI |
416 |
iend = nGroupsInRow |
417 |
jend = nGroupsInCol |
418 |
#else |
419 |
iend = nGroups |
420 |
jend = nGroups |
421 |
#endif |
422 |
|
423 |
!! allocate the groupToGtype and gtypeMaxCutoff here. |
424 |
if(.not.allocated(groupToGtypeRow)) then |
425 |
! allocate(groupToGtype(iend)) |
426 |
allocate(groupToGtypeRow(iend)) |
427 |
else |
428 |
deallocate(groupToGtypeRow) |
429 |
allocate(groupToGtypeRow(iend)) |
430 |
endif |
431 |
if(.not.allocated(groupMaxCutoffRow)) then |
432 |
allocate(groupMaxCutoffRow(iend)) |
433 |
else |
434 |
deallocate(groupMaxCutoffRow) |
435 |
allocate(groupMaxCutoffRow(iend)) |
436 |
end if |
437 |
|
438 |
if(.not.allocated(gtypeMaxCutoffRow)) then |
439 |
allocate(gtypeMaxCutoffRow(iend)) |
440 |
else |
441 |
deallocate(gtypeMaxCutoffRow) |
442 |
allocate(gtypeMaxCutoffRow(iend)) |
443 |
endif |
444 |
|
445 |
|
446 |
#ifdef IS_MPI |
447 |
! We only allocate new storage if we are in MPI because Ncol /= Nrow |
448 |
if(.not.associated(groupToGtypeCol)) then |
449 |
allocate(groupToGtypeCol(jend)) |
450 |
else |
451 |
deallocate(groupToGtypeCol) |
452 |
allocate(groupToGtypeCol(jend)) |
453 |
end if |
454 |
|
455 |
if(.not.associated(groupMaxCutoffCol)) then |
456 |
allocate(groupMaxCutoffCol(jend)) |
457 |
else |
458 |
deallocate(groupMaxCutoffCol) |
459 |
allocate(groupMaxCutoffCol(jend)) |
460 |
end if |
461 |
if(.not.associated(gtypeMaxCutoffCol)) then |
462 |
allocate(gtypeMaxCutoffCol(jend)) |
463 |
else |
464 |
deallocate(gtypeMaxCutoffCol) |
465 |
allocate(gtypeMaxCutoffCol(jend)) |
466 |
end if |
467 |
|
468 |
groupMaxCutoffCol = 0.0_dp |
469 |
gtypeMaxCutoffCol = 0.0_dp |
470 |
|
471 |
#endif |
472 |
groupMaxCutoffRow = 0.0_dp |
473 |
gtypeMaxCutoffRow = 0.0_dp |
474 |
|
475 |
|
476 |
!! first we do a single loop over the cutoff groups to find the |
477 |
!! largest cutoff for any atypes present in this group. We also |
478 |
!! create gtypes at this point. |
479 |
|
480 |
tol = 1.0e-6_dp |
481 |
nGroupTypesRow = 0 |
482 |
nGroupTypesCol = 0 |
483 |
do i = istart, iend |
484 |
n_in_i = groupStartRow(i+1) - groupStartRow(i) |
485 |
groupMaxCutoffRow(i) = 0.0_dp |
486 |
do ia = groupStartRow(i), groupStartRow(i+1)-1 |
487 |
atom1 = groupListRow(ia) |
488 |
#ifdef IS_MPI |
489 |
me_i = atid_row(atom1) |
490 |
#else |
491 |
me_i = atid(atom1) |
492 |
#endif |
493 |
if (atypeMaxCutoff(me_i).gt.groupMaxCutoffRow(i)) then |
494 |
groupMaxCutoffRow(i)=atypeMaxCutoff(me_i) |
495 |
endif |
496 |
enddo |
497 |
if (nGroupTypesRow.eq.0) then |
498 |
nGroupTypesRow = nGroupTypesRow + 1 |
499 |
gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i) |
500 |
groupToGtypeRow(i) = nGroupTypesRow |
501 |
else |
502 |
GtypeFound = .false. |
503 |
do g = 1, nGroupTypesRow |
504 |
if ( abs(groupMaxCutoffRow(i) - gtypeMaxCutoffRow(g)).lt.tol) then |
505 |
groupToGtypeRow(i) = g |
506 |
GtypeFound = .true. |
507 |
endif |
508 |
enddo |
509 |
if (.not.GtypeFound) then |
510 |
nGroupTypesRow = nGroupTypesRow + 1 |
511 |
gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i) |
512 |
groupToGtypeRow(i) = nGroupTypesRow |
513 |
endif |
514 |
endif |
515 |
enddo |
516 |
|
517 |
#ifdef IS_MPI |
518 |
do j = jstart, jend |
519 |
n_in_j = groupStartCol(j+1) - groupStartCol(j) |
520 |
groupMaxCutoffCol(j) = 0.0_dp |
521 |
do ja = groupStartCol(j), groupStartCol(j+1)-1 |
522 |
atom1 = groupListCol(ja) |
523 |
|
524 |
me_j = atid_col(atom1) |
525 |
|
526 |
if (atypeMaxCutoff(me_j).gt.groupMaxCutoffCol(j)) then |
527 |
groupMaxCutoffCol(j)=atypeMaxCutoff(me_j) |
528 |
endif |
529 |
enddo |
530 |
|
531 |
if (nGroupTypesCol.eq.0) then |
532 |
nGroupTypesCol = nGroupTypesCol + 1 |
533 |
gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j) |
534 |
groupToGtypeCol(j) = nGroupTypesCol |
535 |
else |
536 |
GtypeFound = .false. |
537 |
do g = 1, nGroupTypesCol |
538 |
if ( abs(groupMaxCutoffCol(j) - gtypeMaxCutoffCol(g)).lt.tol) then |
539 |
groupToGtypeCol(j) = g |
540 |
GtypeFound = .true. |
541 |
endif |
542 |
enddo |
543 |
if (.not.GtypeFound) then |
544 |
nGroupTypesCol = nGroupTypesCol + 1 |
545 |
gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j) |
546 |
groupToGtypeCol(j) = nGroupTypesCol |
547 |
endif |
548 |
endif |
549 |
enddo |
550 |
|
551 |
#else |
552 |
! Set pointers to information we just found |
553 |
nGroupTypesCol = nGroupTypesRow |
554 |
groupToGtypeCol => groupToGtypeRow |
555 |
gtypeMaxCutoffCol => gtypeMaxCutoffRow |
556 |
groupMaxCutoffCol => groupMaxCutoffRow |
557 |
#endif |
558 |
|
559 |
!! allocate the gtypeCutoffMap here. |
560 |
allocate(gtypeCutoffMap(nGroupTypesRow,nGroupTypesCol)) |
561 |
!! then we do a double loop over all the group TYPES to find the cutoff |
562 |
!! map between groups of two types |
563 |
tradRcut = max(maxval(gtypeMaxCutoffRow),maxval(gtypeMaxCutoffCol)) |
564 |
|
565 |
do i = 1, nGroupTypesRow |
566 |
do j = 1, nGroupTypesCol |
567 |
|
568 |
select case(cutoffPolicy) |
569 |
case(TRADITIONAL_CUTOFF_POLICY) |
570 |
thisRcut = tradRcut |
571 |
case(MIX_CUTOFF_POLICY) |
572 |
thisRcut = 0.5_dp * (gtypeMaxCutoffRow(i) + gtypeMaxCutoffCol(j)) |
573 |
case(MAX_CUTOFF_POLICY) |
574 |
thisRcut = max(gtypeMaxCutoffRow(i), gtypeMaxCutoffCol(j)) |
575 |
case default |
576 |
call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy") |
577 |
return |
578 |
end select |
579 |
gtypeCutoffMap(i,j)%rcut = thisRcut |
580 |
|
581 |
if (thisRcut.gt.largestRcut) largestRcut = thisRcut |
582 |
|
583 |
gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut |
584 |
|
585 |
if (.not.haveSkinThickness) then |
586 |
skinThickness = 1.0_dp |
587 |
endif |
588 |
|
589 |
gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skinThickness)**2 |
590 |
|
591 |
! sanity check |
592 |
|
593 |
if (haveDefaultCutoffs) then |
594 |
if (abs(gtypeCutoffMap(i,j)%rcut - defaultRcut).gt.0.0001) then |
595 |
call handleError("createGtypeCutoffMap", "user-specified rCut does not match computed group Cutoff") |
596 |
endif |
597 |
endif |
598 |
enddo |
599 |
enddo |
600 |
|
601 |
if(allocated(gtypeMaxCutoffRow)) deallocate(gtypeMaxCutoffRow) |
602 |
if(allocated(groupMaxCutoffRow)) deallocate(groupMaxCutoffRow) |
603 |
if(allocated(atypeMaxCutoff)) deallocate(atypeMaxCutoff) |
604 |
#ifdef IS_MPI |
605 |
if(associated(groupMaxCutoffCol)) deallocate(groupMaxCutoffCol) |
606 |
if(associated(gtypeMaxCutoffCol)) deallocate(gtypeMaxCutoffCol) |
607 |
#endif |
608 |
groupMaxCutoffCol => null() |
609 |
gtypeMaxCutoffCol => null() |
610 |
|
611 |
haveGtypeCutoffMap = .true. |
612 |
end subroutine createGtypeCutoffMap |
613 |
|
614 |
subroutine setCutoffs(defRcut, defRsw, defSP, defSF) |
615 |
|
616 |
real(kind=dp),intent(in) :: defRcut, defRsw |
617 |
integer, intent(in) :: defSP, defSF |
618 |
character(len = statusMsgSize) :: errMsg |
619 |
integer :: localError |
620 |
|
621 |
defaultRcut = defRcut |
622 |
defaultRsw = defRsw |
623 |
|
624 |
if (defSP .ne. 0) then |
625 |
defaultDoShiftPot = .true. |
626 |
else |
627 |
defaultDoShiftPot = .false. |
628 |
endif |
629 |
if (defSF .ne. 0) then |
630 |
defaultDoShiftFrc = .true. |
631 |
else |
632 |
defaultDoShiftFrc = .false. |
633 |
endif |
634 |
|
635 |
if (abs(defaultRcut-defaultRsw) .lt. 0.0001) then |
636 |
if (defaultDoShiftFrc) then |
637 |
write(errMsg, *) & |
638 |
'cutoffRadius and switchingRadius are set to the', newline & |
639 |
// tab, 'same value. OpenMD will use shifted force', newline & |
640 |
// tab, 'potentials instead of switching functions.' |
641 |
|
642 |
call handleInfo("setCutoffs", errMsg) |
643 |
else |
644 |
write(errMsg, *) & |
645 |
'cutoffRadius and switchingRadius are set to the', newline & |
646 |
// tab, 'same value. OpenMD will use shifted', newline & |
647 |
// tab, 'potentials instead of switching functions.' |
648 |
|
649 |
call handleInfo("setCutoffs", errMsg) |
650 |
|
651 |
defaultDoShiftPot = .true. |
652 |
endif |
653 |
|
654 |
endif |
655 |
|
656 |
localError = 0 |
657 |
call setLJDefaultCutoff( defaultRcut, defaultDoShiftPot, & |
658 |
defaultDoShiftFrc ) |
659 |
call setElectrostaticCutoffRadius( defaultRcut, defaultRsw ) |
660 |
call setCutoffEAM( defaultRcut ) |
661 |
call setCutoffSC( defaultRcut ) |
662 |
call setMnMDefaultCutoff( defaultRcut, defaultDoShiftPot, & |
663 |
defaultDoShiftFrc ) |
664 |
call set_switch(defaultRsw, defaultRcut) |
665 |
call setHmatDangerousRcutValue(defaultRcut) |
666 |
|
667 |
haveDefaultCutoffs = .true. |
668 |
haveGtypeCutoffMap = .false. |
669 |
|
670 |
end subroutine setCutoffs |
671 |
|
672 |
subroutine cWasLame() |
673 |
|
674 |
VisitCutoffsAfterComputing = .true. |
675 |
return |
676 |
|
677 |
end subroutine cWasLame |
678 |
|
679 |
subroutine setCutoffPolicy(cutPolicy) |
680 |
|
681 |
integer, intent(in) :: cutPolicy |
682 |
|
683 |
cutoffPolicy = cutPolicy |
684 |
haveCutoffPolicy = .true. |
685 |
haveGtypeCutoffMap = .false. |
686 |
|
687 |
end subroutine setCutoffPolicy |
688 |
|
689 |
subroutine setBoxDipole() |
690 |
|
691 |
do_box_dipole = .true. |
692 |
|
693 |
end subroutine setBoxDipole |
694 |
|
695 |
subroutine getBoxDipole( box_dipole ) |
696 |
|
697 |
real(kind=dp), intent(inout), dimension(3) :: box_dipole |
698 |
|
699 |
box_dipole = boxDipole |
700 |
|
701 |
end subroutine getBoxDipole |
702 |
|
703 |
subroutine setElectrostaticMethod( thisESM ) |
704 |
|
705 |
integer, intent(in) :: thisESM |
706 |
|
707 |
electrostaticSummationMethod = thisESM |
708 |
haveElectrostaticSummationMethod = .true. |
709 |
|
710 |
end subroutine setElectrostaticMethod |
711 |
|
712 |
subroutine setSkinThickness( thisSkin ) |
713 |
|
714 |
real(kind=dp), intent(in) :: thisSkin |
715 |
|
716 |
skinThickness = thisSkin |
717 |
haveSkinThickness = .true. |
718 |
haveGtypeCutoffMap = .false. |
719 |
|
720 |
end subroutine setSkinThickness |
721 |
|
722 |
subroutine setSimVariables() |
723 |
SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms() |
724 |
SIM_uses_EAM = SimUsesEAM() |
725 |
SIM_requires_postpair_calc = SimRequiresPostpairCalc() |
726 |
SIM_requires_prepair_calc = SimRequiresPrepairCalc() |
727 |
SIM_uses_PBC = SimUsesPBC() |
728 |
SIM_uses_SC = SimUsesSC() |
729 |
SIM_uses_AtomicVirial = SimUsesAtomicVirial() |
730 |
|
731 |
haveSIMvariables = .true. |
732 |
|
733 |
return |
734 |
end subroutine setSimVariables |
735 |
|
736 |
subroutine doReadyCheck(error) |
737 |
integer, intent(out) :: error |
738 |
integer :: myStatus |
739 |
|
740 |
error = 0 |
741 |
|
742 |
if (.not. haveInteractionHash) then |
743 |
call createInteractionHash() |
744 |
endif |
745 |
|
746 |
if (.not. haveGtypeCutoffMap) then |
747 |
call createGtypeCutoffMap() |
748 |
endif |
749 |
|
750 |
if (VisitCutoffsAfterComputing) then |
751 |
call set_switch(largestRcut, largestRcut) |
752 |
call setHmatDangerousRcutValue(largestRcut) |
753 |
call setCutoffEAM(largestRcut) |
754 |
call setCutoffSC(largestRcut) |
755 |
VisitCutoffsAfterComputing = .false. |
756 |
endif |
757 |
|
758 |
if (.not. haveSIMvariables) then |
759 |
call setSimVariables() |
760 |
endif |
761 |
|
762 |
if (.not. haveNeighborList) then |
763 |
write(default_error, *) 'neighbor list has not been initialized in doForces!' |
764 |
error = -1 |
765 |
return |
766 |
end if |
767 |
|
768 |
if (.not. haveSaneForceField) then |
769 |
write(default_error, *) 'Force Field is not sane in doForces!' |
770 |
error = -1 |
771 |
return |
772 |
end if |
773 |
|
774 |
#ifdef IS_MPI |
775 |
if (.not. isMPISimSet()) then |
776 |
write(default_error,*) "ERROR: mpiSimulation has not been initialized!" |
777 |
error = -1 |
778 |
return |
779 |
endif |
780 |
#endif |
781 |
return |
782 |
end subroutine doReadyCheck |
783 |
|
784 |
|
785 |
subroutine init_FF(thisStat) |
786 |
|
787 |
integer, intent(out) :: thisStat |
788 |
integer :: my_status, nMatches |
789 |
integer, pointer :: MatchList(:) => null() |
790 |
|
791 |
!! assume things are copacetic, unless they aren't |
792 |
thisStat = 0 |
793 |
|
794 |
!! init_FF is called *after* all of the atom types have been |
795 |
!! defined in atype_module using the new_atype subroutine. |
796 |
!! |
797 |
!! this will scan through the known atypes and figure out what |
798 |
!! interactions are used by the force field. |
799 |
|
800 |
FF_uses_DirectionalAtoms = .false. |
801 |
FF_uses_Dipoles = .false. |
802 |
FF_uses_GayBerne = .false. |
803 |
FF_uses_EAM = .false. |
804 |
FF_uses_SC = .false. |
805 |
|
806 |
call getMatchingElementList(atypes, "is_Directional", .true., & |
807 |
nMatches, MatchList) |
808 |
if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true. |
809 |
|
810 |
call getMatchingElementList(atypes, "is_Dipole", .true., & |
811 |
nMatches, MatchList) |
812 |
if (nMatches .gt. 0) FF_uses_Dipoles = .true. |
813 |
|
814 |
call getMatchingElementList(atypes, "is_GayBerne", .true., & |
815 |
nMatches, MatchList) |
816 |
if (nMatches .gt. 0) FF_uses_GayBerne = .true. |
817 |
|
818 |
call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList) |
819 |
if (nMatches .gt. 0) FF_uses_EAM = .true. |
820 |
|
821 |
call getMatchingElementList(atypes, "is_SC", .true., nMatches, MatchList) |
822 |
if (nMatches .gt. 0) FF_uses_SC = .true. |
823 |
|
824 |
|
825 |
haveSaneForceField = .true. |
826 |
|
827 |
|
828 |
if (.not. haveNeighborList) then |
829 |
!! Create neighbor lists |
830 |
call expandNeighborList(nLocal, my_status) |
831 |
if (my_Status /= 0) then |
832 |
write(default_error,*) "SimSetup: ExpandNeighborList returned error." |
833 |
thisStat = -1 |
834 |
return |
835 |
endif |
836 |
haveNeighborList = .true. |
837 |
endif |
838 |
|
839 |
!!$ ! get the C-side pointers |
840 |
!!$ gsfunptr = getSigmaFunc() |
841 |
!!$ gefunptr = getEpsilonFunc() |
842 |
!!$ dljpsubptr = getLJPfunc() |
843 |
!!$ ! attach to fortran |
844 |
!!$ call c_f_procpointer(gsfunptr, getSigma) |
845 |
!!$ call c_f_procpointer(gefunptr, getEpsilon) |
846 |
!!$ call c_f_procpointer(dljpsubptr, do_lj_pair) |
847 |
|
848 |
end subroutine init_FF |
849 |
|
850 |
|
851 |
!! Does force loop over i,j pairs. Calls do_pair to calculates forces. |
852 |
!-------------------------------------------------------------> |
853 |
subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, particle_pot, & |
854 |
error) |
855 |
!! Position array provided by C, dimensioned by getNlocal |
856 |
real ( kind = dp ), dimension(3, nLocal) :: q |
857 |
!! molecular center-of-mass position array |
858 |
real ( kind = dp ), dimension(3, nGroups) :: q_group |
859 |
!! Rotation Matrix for each long range particle in simulation. |
860 |
real( kind = dp), dimension(9, nLocal) :: A |
861 |
!! Unit vectors for dipoles (lab frame) |
862 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
863 |
!! Force array provided by C, dimensioned by getNlocal |
864 |
real ( kind = dp ), dimension(3,nLocal) :: f |
865 |
!! Torsion array provided by C, dimensioned by getNlocal |
866 |
real( kind = dp ), dimension(3,nLocal) :: t |
867 |
|
868 |
!! Stress Tensor |
869 |
real( kind = dp), dimension(9) :: tau |
870 |
real ( kind = dp ),dimension(LR_POT_TYPES) :: pot |
871 |
real( kind = dp ), dimension(nLocal) :: particle_pot |
872 |
|
873 |
logical :: in_switching_region |
874 |
#ifdef IS_MPI |
875 |
real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local |
876 |
integer :: nAtomsInRow |
877 |
integer :: nAtomsInCol |
878 |
integer :: nprocs |
879 |
integer :: nGroupsInRow |
880 |
integer :: nGroupsInCol |
881 |
#endif |
882 |
integer :: natoms |
883 |
logical :: update_nlist |
884 |
integer :: i, j, jstart, jend, jnab |
885 |
integer :: istart, iend |
886 |
integer :: ia, jb, atom1, atom2 |
887 |
integer :: nlist |
888 |
real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij |
889 |
real( kind = DP ) :: sw, dswdr, swderiv, mf |
890 |
real( kind = DP ) :: rVal |
891 |
real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag |
892 |
real(kind=dp) :: rfpot, mu_i |
893 |
real(kind=dp):: rCut |
894 |
integer :: me_i, me_j, n_in_i, n_in_j, iG, j1 |
895 |
logical :: is_dp_i |
896 |
integer :: neighborListSize |
897 |
integer :: listerror, error |
898 |
integer :: localError |
899 |
integer :: propPack_i, propPack_j |
900 |
integer :: loopStart, loopEnd, loop |
901 |
integer :: iHash, jHash |
902 |
integer :: i1, topoDist |
903 |
|
904 |
!! the variables for the box dipole moment |
905 |
#ifdef IS_MPI |
906 |
integer :: pChgCount_local |
907 |
integer :: nChgCount_local |
908 |
real(kind=dp) :: pChg_local |
909 |
real(kind=dp) :: nChg_local |
910 |
real(kind=dp), dimension(3) :: pChgPos_local |
911 |
real(kind=dp), dimension(3) :: nChgPos_local |
912 |
real(kind=dp), dimension(3) :: dipVec_local |
913 |
#endif |
914 |
integer :: pChgCount |
915 |
integer :: nChgCount |
916 |
real(kind=dp) :: pChg |
917 |
real(kind=dp) :: nChg |
918 |
real(kind=dp) :: chg_value |
919 |
real(kind=dp), dimension(3) :: pChgPos |
920 |
real(kind=dp), dimension(3) :: nChgPos |
921 |
real(kind=dp), dimension(3) :: dipVec |
922 |
real(kind=dp), dimension(3) :: chgVec |
923 |
real(kind=dp) :: skch |
924 |
|
925 |
!! initialize box dipole variables |
926 |
if (do_box_dipole) then |
927 |
#ifdef IS_MPI |
928 |
pChg_local = 0.0_dp |
929 |
nChg_local = 0.0_dp |
930 |
pChgCount_local = 0 |
931 |
nChgCount_local = 0 |
932 |
do i=1, 3 |
933 |
pChgPos_local = 0.0_dp |
934 |
nChgPos_local = 0.0_dp |
935 |
dipVec_local = 0.0_dp |
936 |
enddo |
937 |
#endif |
938 |
pChg = 0.0_dp |
939 |
nChg = 0.0_dp |
940 |
pChgCount = 0 |
941 |
nChgCount = 0 |
942 |
chg_value = 0.0_dp |
943 |
|
944 |
do i=1, 3 |
945 |
pChgPos(i) = 0.0_dp |
946 |
nChgPos(i) = 0.0_dp |
947 |
dipVec(i) = 0.0_dp |
948 |
chgVec(i) = 0.0_dp |
949 |
boxDipole(i) = 0.0_dp |
950 |
enddo |
951 |
endif |
952 |
|
953 |
!! initialize local variables |
954 |
|
955 |
#ifdef IS_MPI |
956 |
pot_local = 0.0_dp |
957 |
nAtomsInRow = getNatomsInRow(plan_atom_row) |
958 |
nAtomsInCol = getNatomsInCol(plan_atom_col) |
959 |
nGroupsInRow = getNgroupsInRow(plan_group_row) |
960 |
nGroupsInCol = getNgroupsInCol(plan_group_col) |
961 |
#else |
962 |
natoms = nlocal |
963 |
#endif |
964 |
|
965 |
call doReadyCheck(localError) |
966 |
if ( localError .ne. 0 ) then |
967 |
call handleError("do_force_loop", "Not Initialized") |
968 |
error = -1 |
969 |
return |
970 |
end if |
971 |
call zero_work_arrays() |
972 |
|
973 |
! Gather all information needed by all force loops: |
974 |
|
975 |
#ifdef IS_MPI |
976 |
|
977 |
call gather(q, q_Row, plan_atom_row_3d) |
978 |
call gather(q, q_Col, plan_atom_col_3d) |
979 |
|
980 |
call gather(q_group, q_group_Row, plan_group_row_3d) |
981 |
call gather(q_group, q_group_Col, plan_group_col_3d) |
982 |
|
983 |
if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then |
984 |
call gather(eFrame, eFrame_Row, plan_atom_row_rotation) |
985 |
call gather(eFrame, eFrame_Col, plan_atom_col_rotation) |
986 |
|
987 |
call gather(A, A_Row, plan_atom_row_rotation) |
988 |
call gather(A, A_Col, plan_atom_col_rotation) |
989 |
endif |
990 |
|
991 |
#endif |
992 |
|
993 |
!! Begin force loop timing: |
994 |
#ifdef PROFILE |
995 |
call cpu_time(forceTimeInitial) |
996 |
nloops = nloops + 1 |
997 |
#endif |
998 |
|
999 |
loopEnd = PAIR_LOOP |
1000 |
if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then |
1001 |
loopStart = PREPAIR_LOOP |
1002 |
else |
1003 |
loopStart = PAIR_LOOP |
1004 |
endif |
1005 |
|
1006 |
do loop = loopStart, loopEnd |
1007 |
|
1008 |
! See if we need to update neighbor lists |
1009 |
! (but only on the first time through): |
1010 |
if (loop .eq. loopStart) then |
1011 |
#ifdef IS_MPI |
1012 |
call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, & |
1013 |
update_nlist) |
1014 |
#else |
1015 |
call checkNeighborList(nGroups, q_group, skinThickness, & |
1016 |
update_nlist) |
1017 |
#endif |
1018 |
endif |
1019 |
|
1020 |
if (update_nlist) then |
1021 |
!! save current configuration and construct neighbor list |
1022 |
#ifdef IS_MPI |
1023 |
call saveNeighborList(nGroupsInRow, q_group_row) |
1024 |
#else |
1025 |
call saveNeighborList(nGroups, q_group) |
1026 |
#endif |
1027 |
neighborListSize = size(list) |
1028 |
nlist = 0 |
1029 |
endif |
1030 |
|
1031 |
istart = 1 |
1032 |
#ifdef IS_MPI |
1033 |
iend = nGroupsInRow |
1034 |
#else |
1035 |
iend = nGroups - 1 |
1036 |
#endif |
1037 |
outer: do i = istart, iend |
1038 |
|
1039 |
if (update_nlist) point(i) = nlist + 1 |
1040 |
|
1041 |
n_in_i = groupStartRow(i+1) - groupStartRow(i) |
1042 |
|
1043 |
if (update_nlist) then |
1044 |
#ifdef IS_MPI |
1045 |
jstart = 1 |
1046 |
jend = nGroupsInCol |
1047 |
#else |
1048 |
jstart = i+1 |
1049 |
jend = nGroups |
1050 |
#endif |
1051 |
else |
1052 |
jstart = point(i) |
1053 |
jend = point(i+1) - 1 |
1054 |
! make sure group i has neighbors |
1055 |
if (jstart .gt. jend) cycle outer |
1056 |
endif |
1057 |
|
1058 |
do jnab = jstart, jend |
1059 |
if (update_nlist) then |
1060 |
j = jnab |
1061 |
else |
1062 |
j = list(jnab) |
1063 |
endif |
1064 |
|
1065 |
#ifdef IS_MPI |
1066 |
me_j = atid_col(j) |
1067 |
call get_interatomic_vector(q_group_Row(:,i), & |
1068 |
q_group_Col(:,j), d_grp, rgrpsq) |
1069 |
#else |
1070 |
me_j = atid(j) |
1071 |
call get_interatomic_vector(q_group(:,i), & |
1072 |
q_group(:,j), d_grp, rgrpsq) |
1073 |
#endif |
1074 |
|
1075 |
if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then |
1076 |
if (update_nlist) then |
1077 |
nlist = nlist + 1 |
1078 |
|
1079 |
if (nlist > neighborListSize) then |
1080 |
#ifdef IS_MPI |
1081 |
call expandNeighborList(nGroupsInRow, listerror) |
1082 |
#else |
1083 |
call expandNeighborList(nGroups, listerror) |
1084 |
#endif |
1085 |
if (listerror /= 0) then |
1086 |
error = -1 |
1087 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
1088 |
return |
1089 |
end if |
1090 |
neighborListSize = size(list) |
1091 |
endif |
1092 |
|
1093 |
list(nlist) = j |
1094 |
endif |
1095 |
|
1096 |
if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then |
1097 |
|
1098 |
rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut |
1099 |
if (loop .eq. PAIR_LOOP) then |
1100 |
vij = 0.0_dp |
1101 |
fij(1) = 0.0_dp |
1102 |
fij(2) = 0.0_dp |
1103 |
fij(3) = 0.0_dp |
1104 |
endif |
1105 |
|
1106 |
call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region) |
1107 |
|
1108 |
n_in_j = groupStartCol(j+1) - groupStartCol(j) |
1109 |
|
1110 |
do ia = groupStartRow(i), groupStartRow(i+1)-1 |
1111 |
|
1112 |
atom1 = groupListRow(ia) |
1113 |
|
1114 |
inner: do jb = groupStartCol(j), groupStartCol(j+1)-1 |
1115 |
|
1116 |
atom2 = groupListCol(jb) |
1117 |
|
1118 |
if (skipThisPair(atom1, atom2)) cycle inner |
1119 |
|
1120 |
if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then |
1121 |
d_atm(1) = d_grp(1) |
1122 |
d_atm(2) = d_grp(2) |
1123 |
d_atm(3) = d_grp(3) |
1124 |
ratmsq = rgrpsq |
1125 |
else |
1126 |
#ifdef IS_MPI |
1127 |
call get_interatomic_vector(q_Row(:,atom1), & |
1128 |
q_Col(:,atom2), d_atm, ratmsq) |
1129 |
#else |
1130 |
call get_interatomic_vector(q(:,atom1), & |
1131 |
q(:,atom2), d_atm, ratmsq) |
1132 |
#endif |
1133 |
endif |
1134 |
|
1135 |
topoDist = getTopoDistance(atom1, atom2) |
1136 |
|
1137 |
if (loop .eq. PREPAIR_LOOP) then |
1138 |
#ifdef IS_MPI |
1139 |
call do_prepair(atom1, atom2, ratmsq, d_atm, sw, & |
1140 |
rgrpsq, d_grp, rCut, & |
1141 |
eFrame, A, f, t, pot_local) |
1142 |
#else |
1143 |
call do_prepair(atom1, atom2, ratmsq, d_atm, sw, & |
1144 |
rgrpsq, d_grp, rCut, & |
1145 |
eFrame, A, f, t, pot) |
1146 |
#endif |
1147 |
else |
1148 |
#ifdef IS_MPI |
1149 |
call do_pair(atom1, atom2, ratmsq, d_atm, sw, & |
1150 |
eFrame, A, f, t, pot_local, particle_pot, vpair, & |
1151 |
fpair, d_grp, rgrp, rCut, topoDist) |
1152 |
! particle_pot will be accumulated from row & column |
1153 |
! arrays later |
1154 |
#else |
1155 |
call do_pair(atom1, atom2, ratmsq, d_atm, sw, & |
1156 |
eFrame, A, f, t, pot, particle_pot, vpair, & |
1157 |
fpair, d_grp, rgrp, rCut, topoDist) |
1158 |
#endif |
1159 |
vij = vij + vpair |
1160 |
fij(1) = fij(1) + fpair(1) |
1161 |
fij(2) = fij(2) + fpair(2) |
1162 |
fij(3) = fij(3) + fpair(3) |
1163 |
call add_stress_tensor(d_atm, fpair, tau) |
1164 |
endif |
1165 |
enddo inner |
1166 |
enddo |
1167 |
|
1168 |
if (loop .eq. PAIR_LOOP) then |
1169 |
if (in_switching_region) then |
1170 |
swderiv = vij*dswdr/rgrp |
1171 |
fg = swderiv*d_grp |
1172 |
|
1173 |
fij(1) = fij(1) + fg(1) |
1174 |
fij(2) = fij(2) + fg(2) |
1175 |
fij(3) = fij(3) + fg(3) |
1176 |
|
1177 |
if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then |
1178 |
call add_stress_tensor(d_atm, fg, tau) |
1179 |
endif |
1180 |
|
1181 |
do ia=groupStartRow(i), groupStartRow(i+1)-1 |
1182 |
atom1=groupListRow(ia) |
1183 |
mf = mfactRow(atom1) |
1184 |
! fg is the force on atom ia due to cutoff group's |
1185 |
! presence in switching region |
1186 |
fg = swderiv*d_grp*mf |
1187 |
#ifdef IS_MPI |
1188 |
f_Row(1,atom1) = f_Row(1,atom1) + fg(1) |
1189 |
f_Row(2,atom1) = f_Row(2,atom1) + fg(2) |
1190 |
f_Row(3,atom1) = f_Row(3,atom1) + fg(3) |
1191 |
#else |
1192 |
f(1,atom1) = f(1,atom1) + fg(1) |
1193 |
f(2,atom1) = f(2,atom1) + fg(2) |
1194 |
f(3,atom1) = f(3,atom1) + fg(3) |
1195 |
#endif |
1196 |
if (n_in_i .gt. 1) then |
1197 |
if (SIM_uses_AtomicVirial) then |
1198 |
! find the distance between the atom |
1199 |
! and the center of the cutoff group: |
1200 |
#ifdef IS_MPI |
1201 |
call get_interatomic_vector(q_Row(:,atom1), & |
1202 |
q_group_Row(:,i), dag, rag) |
1203 |
#else |
1204 |
call get_interatomic_vector(q(:,atom1), & |
1205 |
q_group(:,i), dag, rag) |
1206 |
#endif |
1207 |
call add_stress_tensor(dag,fg,tau) |
1208 |
endif |
1209 |
endif |
1210 |
enddo |
1211 |
|
1212 |
do jb=groupStartCol(j), groupStartCol(j+1)-1 |
1213 |
atom2=groupListCol(jb) |
1214 |
mf = mfactCol(atom2) |
1215 |
! fg is the force on atom jb due to cutoff group's |
1216 |
! presence in switching region |
1217 |
fg = -swderiv*d_grp*mf |
1218 |
#ifdef IS_MPI |
1219 |
f_Col(1,atom2) = f_Col(1,atom2) + fg(1) |
1220 |
f_Col(2,atom2) = f_Col(2,atom2) + fg(2) |
1221 |
f_Col(3,atom2) = f_Col(3,atom2) + fg(3) |
1222 |
#else |
1223 |
f(1,atom2) = f(1,atom2) + fg(1) |
1224 |
f(2,atom2) = f(2,atom2) + fg(2) |
1225 |
f(3,atom2) = f(3,atom2) + fg(3) |
1226 |
#endif |
1227 |
if (n_in_j .gt. 1) then |
1228 |
if (SIM_uses_AtomicVirial) then |
1229 |
! find the distance between the atom |
1230 |
! and the center of the cutoff group: |
1231 |
#ifdef IS_MPI |
1232 |
call get_interatomic_vector(q_Col(:,atom2), & |
1233 |
q_group_Col(:,j), dag, rag) |
1234 |
#else |
1235 |
call get_interatomic_vector(q(:,atom2), & |
1236 |
q_group(:,j), dag, rag) |
1237 |
#endif |
1238 |
call add_stress_tensor(dag,fg,tau) |
1239 |
endif |
1240 |
endif |
1241 |
enddo |
1242 |
endif |
1243 |
!if (.not.SIM_uses_AtomicVirial) then |
1244 |
! call add_stress_tensor(d_grp, fij, tau) |
1245 |
!endif |
1246 |
endif |
1247 |
endif |
1248 |
endif |
1249 |
enddo |
1250 |
|
1251 |
enddo outer |
1252 |
|
1253 |
if (update_nlist) then |
1254 |
#ifdef IS_MPI |
1255 |
point(nGroupsInRow + 1) = nlist + 1 |
1256 |
#else |
1257 |
point(nGroups) = nlist + 1 |
1258 |
#endif |
1259 |
if (loop .eq. PREPAIR_LOOP) then |
1260 |
! we just did the neighbor list update on the first |
1261 |
! pass, so we don't need to do it |
1262 |
! again on the second pass |
1263 |
update_nlist = .false. |
1264 |
endif |
1265 |
endif |
1266 |
|
1267 |
if (loop .eq. PREPAIR_LOOP) then |
1268 |
#ifdef IS_MPI |
1269 |
call do_preforce(nlocal, pot_local, particle_pot) |
1270 |
#else |
1271 |
call do_preforce(nlocal, pot, particle_pot) |
1272 |
#endif |
1273 |
endif |
1274 |
|
1275 |
enddo |
1276 |
|
1277 |
!! Do timing |
1278 |
#ifdef PROFILE |
1279 |
call cpu_time(forceTimeFinal) |
1280 |
forceTime = forceTime + forceTimeFinal - forceTimeInitial |
1281 |
#endif |
1282 |
|
1283 |
#ifdef IS_MPI |
1284 |
!!distribute forces |
1285 |
|
1286 |
f_temp = 0.0_dp |
1287 |
call scatter(f_Row,f_temp,plan_atom_row_3d) |
1288 |
do i = 1,nlocal |
1289 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
1290 |
end do |
1291 |
|
1292 |
f_temp = 0.0_dp |
1293 |
call scatter(f_Col,f_temp,plan_atom_col_3d) |
1294 |
do i = 1,nlocal |
1295 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
1296 |
end do |
1297 |
|
1298 |
if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then |
1299 |
t_temp = 0.0_dp |
1300 |
call scatter(t_Row,t_temp,plan_atom_row_3d) |
1301 |
do i = 1,nlocal |
1302 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
1303 |
end do |
1304 |
t_temp = 0.0_dp |
1305 |
call scatter(t_Col,t_temp,plan_atom_col_3d) |
1306 |
|
1307 |
do i = 1,nlocal |
1308 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
1309 |
end do |
1310 |
endif |
1311 |
|
1312 |
! scatter/gather pot_row into the members of my column |
1313 |
do i = 1,LR_POT_TYPES |
1314 |
call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row) |
1315 |
end do |
1316 |
! scatter/gather pot_local into all other procs |
1317 |
! add resultant to get total pot |
1318 |
do i = 1, nlocal |
1319 |
pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) & |
1320 |
+ pot_Temp(1:LR_POT_TYPES,i) |
1321 |
enddo |
1322 |
|
1323 |
do i = 1,LR_POT_TYPES |
1324 |
particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal) |
1325 |
enddo |
1326 |
|
1327 |
pot_Temp = 0.0_DP |
1328 |
|
1329 |
do i = 1,LR_POT_TYPES |
1330 |
call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col) |
1331 |
end do |
1332 |
|
1333 |
do i = 1, nlocal |
1334 |
pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)& |
1335 |
+ pot_Temp(1:LR_POT_TYPES,i) |
1336 |
enddo |
1337 |
|
1338 |
do i = 1,LR_POT_TYPES |
1339 |
particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal) |
1340 |
enddo |
1341 |
|
1342 |
ppot_Temp = 0.0_DP |
1343 |
|
1344 |
call scatter(ppot_Row(:), ppot_Temp(:), plan_atom_row) |
1345 |
do i = 1, nlocal |
1346 |
particle_pot(i) = particle_pot(i) + ppot_Temp(i) |
1347 |
enddo |
1348 |
|
1349 |
ppot_Temp = 0.0_DP |
1350 |
|
1351 |
call scatter(ppot_Col(:), ppot_Temp(:), plan_atom_col) |
1352 |
do i = 1, nlocal |
1353 |
particle_pot(i) = particle_pot(i) + ppot_Temp(i) |
1354 |
enddo |
1355 |
|
1356 |
#endif |
1357 |
|
1358 |
if (SIM_requires_postpair_calc) then |
1359 |
do i = 1, nlocal |
1360 |
|
1361 |
! we loop only over the local atoms, so we don't need row and column |
1362 |
! lookups for the types |
1363 |
|
1364 |
me_i = atid(i) |
1365 |
|
1366 |
! is the atom electrostatic? See if it would have an |
1367 |
! electrostatic interaction with itself |
1368 |
iHash = InteractionHash(me_i,me_i) |
1369 |
|
1370 |
if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then |
1371 |
|
1372 |
! loop over the excludes to accumulate charge in the |
1373 |
! cutoff sphere that we've left out of the normal pair loop |
1374 |
skch = 0.0_dp |
1375 |
|
1376 |
do i1 = 1, nSkipsForLocalAtom(i) |
1377 |
j = skipsForLocalAtom(i, i1) |
1378 |
me_j = atid(j) |
1379 |
jHash = InteractionHash(me_i,me_j) |
1380 |
if ( iand(jHash, ELECTROSTATIC_PAIR).ne.0 ) then |
1381 |
skch = skch + getCharge(me_j) |
1382 |
endif |
1383 |
enddo |
1384 |
|
1385 |
#ifdef IS_MPI |
1386 |
call self_self(i, eFrame, skch, pot_local(ELECTROSTATIC_POT), t) |
1387 |
#else |
1388 |
call self_self(i, eFrame, skch, pot(ELECTROSTATIC_POT), t) |
1389 |
#endif |
1390 |
endif |
1391 |
|
1392 |
|
1393 |
if (electrostaticSummationMethod.eq.REACTION_FIELD) then |
1394 |
|
1395 |
! loop over the excludes to accumulate RF stuff we've |
1396 |
! left out of the normal pair loop |
1397 |
|
1398 |
do i1 = 1, nSkipsForLocalAtom(i) |
1399 |
j = skipsForLocalAtom(i, i1) |
1400 |
|
1401 |
! prevent overcounting of the skips |
1402 |
if (i.lt.j) then |
1403 |
call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq) |
1404 |
rVal = sqrt(ratmsq) |
1405 |
call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region) |
1406 |
#ifdef IS_MPI |
1407 |
call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, & |
1408 |
vpair, pot_local(ELECTROSTATIC_POT), f, t) |
1409 |
#else |
1410 |
call rf_self_excludes(i, j, sw, 1.0_dp, eFrame, d_atm, rVal, & |
1411 |
vpair, pot(ELECTROSTATIC_POT), f, t) |
1412 |
#endif |
1413 |
endif |
1414 |
enddo |
1415 |
endif |
1416 |
|
1417 |
if (do_box_dipole) then |
1418 |
#ifdef IS_MPI |
1419 |
call accumulate_box_dipole(i, eFrame, q(:,i), pChg_local, & |
1420 |
nChg_local, pChgPos_local, nChgPos_local, dipVec_local, & |
1421 |
pChgCount_local, nChgCount_local) |
1422 |
#else |
1423 |
call accumulate_box_dipole(i, eFrame, q(:,i), pChg, nChg, & |
1424 |
pChgPos, nChgPos, dipVec, pChgCount, nChgCount) |
1425 |
#endif |
1426 |
endif |
1427 |
enddo |
1428 |
endif |
1429 |
|
1430 |
#ifdef IS_MPI |
1431 |
#ifdef SINGLE_PRECISION |
1432 |
call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, & |
1433 |
mpi_comm_world,mpi_err) |
1434 |
#else |
1435 |
call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, & |
1436 |
mpi_sum, mpi_comm_world,mpi_err) |
1437 |
#endif |
1438 |
|
1439 |
if (do_box_dipole) then |
1440 |
|
1441 |
#ifdef SINGLE_PRECISION |
1442 |
call mpi_allreduce(pChg_local, pChg, 1, mpi_real, mpi_sum, & |
1443 |
mpi_comm_world, mpi_err) |
1444 |
call mpi_allreduce(nChg_local, nChg, 1, mpi_real, mpi_sum, & |
1445 |
mpi_comm_world, mpi_err) |
1446 |
call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer, mpi_sum,& |
1447 |
mpi_comm_world, mpi_err) |
1448 |
call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer, mpi_sum,& |
1449 |
mpi_comm_world, mpi_err) |
1450 |
call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_real, mpi_sum, & |
1451 |
mpi_comm_world, mpi_err) |
1452 |
call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_real, mpi_sum, & |
1453 |
mpi_comm_world, mpi_err) |
1454 |
call mpi_allreduce(dipVec_local, dipVec, 3, mpi_real, mpi_sum, & |
1455 |
mpi_comm_world, mpi_err) |
1456 |
#else |
1457 |
call mpi_allreduce(pChg_local, pChg, 1, mpi_double_precision, mpi_sum, & |
1458 |
mpi_comm_world, mpi_err) |
1459 |
call mpi_allreduce(nChg_local, nChg, 1, mpi_double_precision, mpi_sum, & |
1460 |
mpi_comm_world, mpi_err) |
1461 |
call mpi_allreduce(pChgCount_local, pChgCount, 1, mpi_integer,& |
1462 |
mpi_sum, mpi_comm_world, mpi_err) |
1463 |
call mpi_allreduce(nChgCount_local, nChgCount, 1, mpi_integer,& |
1464 |
mpi_sum, mpi_comm_world, mpi_err) |
1465 |
call mpi_allreduce(pChgPos_local, pChgPos, 3, mpi_double_precision, & |
1466 |
mpi_sum, mpi_comm_world, mpi_err) |
1467 |
call mpi_allreduce(nChgPos_local, nChgPos, 3, mpi_double_precision, & |
1468 |
mpi_sum, mpi_comm_world, mpi_err) |
1469 |
call mpi_allreduce(dipVec_local, dipVec, 3, mpi_double_precision, & |
1470 |
mpi_sum, mpi_comm_world, mpi_err) |
1471 |
#endif |
1472 |
|
1473 |
endif |
1474 |
|
1475 |
#endif |
1476 |
|
1477 |
if (do_box_dipole) then |
1478 |
! first load the accumulated dipole moment (if dipoles were present) |
1479 |
boxDipole(1) = dipVec(1) |
1480 |
boxDipole(2) = dipVec(2) |
1481 |
boxDipole(3) = dipVec(3) |
1482 |
|
1483 |
! now include the dipole moment due to charges |
1484 |
! use the lesser of the positive and negative charge totals |
1485 |
if (nChg .le. pChg) then |
1486 |
chg_value = nChg |
1487 |
else |
1488 |
chg_value = pChg |
1489 |
endif |
1490 |
|
1491 |
! find the average positions |
1492 |
if (pChgCount .gt. 0 .and. nChgCount .gt. 0) then |
1493 |
pChgPos = pChgPos / pChgCount |
1494 |
nChgPos = nChgPos / nChgCount |
1495 |
endif |
1496 |
|
1497 |
! dipole is from the negative to the positive (physics notation) |
1498 |
chgVec(1) = pChgPos(1) - nChgPos(1) |
1499 |
chgVec(2) = pChgPos(2) - nChgPos(2) |
1500 |
chgVec(3) = pChgPos(3) - nChgPos(3) |
1501 |
|
1502 |
boxDipole(1) = boxDipole(1) + chgVec(1) * chg_value |
1503 |
boxDipole(2) = boxDipole(2) + chgVec(2) * chg_value |
1504 |
boxDipole(3) = boxDipole(3) + chgVec(3) * chg_value |
1505 |
|
1506 |
endif |
1507 |
|
1508 |
end subroutine do_force_loop |
1509 |
|
1510 |
subroutine do_pair(i, j, rijsq, d, sw, & |
1511 |
eFrame, A, f, t, pot, particle_pot, vpair, & |
1512 |
fpair, d_grp, r_grp, rCut, topoDist) |
1513 |
|
1514 |
real( kind = dp ) :: vpair, sw |
1515 |
real( kind = dp ), dimension(LR_POT_TYPES) :: pot |
1516 |
real( kind = dp ), dimension(nLocal) :: particle_pot |
1517 |
real( kind = dp ), dimension(3) :: fpair |
1518 |
real( kind = dp ), dimension(nLocal) :: mfact |
1519 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
1520 |
real( kind = dp ), dimension(9,nLocal) :: A |
1521 |
real( kind = dp ), dimension(3,nLocal) :: f |
1522 |
real( kind = dp ), dimension(3,nLocal) :: t |
1523 |
|
1524 |
integer, intent(in) :: i, j |
1525 |
real ( kind = dp ), intent(inout) :: rijsq |
1526 |
real ( kind = dp ), intent(inout) :: r_grp |
1527 |
real ( kind = dp ), intent(inout) :: d(3) |
1528 |
real ( kind = dp ), intent(inout) :: d_grp(3) |
1529 |
real ( kind = dp ), intent(inout) :: rCut |
1530 |
integer, intent(inout) :: topoDist |
1531 |
real ( kind = dp ) :: r, pair_pot, vdwMult, electroMult |
1532 |
real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx |
1533 |
|
1534 |
real( kind = dp), dimension(3) :: f1, t1, t2 |
1535 |
real( kind = dp), dimension(9) :: A1, A2, eF1, eF2 |
1536 |
real( kind = dp) :: dfrhodrho_i, dfrhodrho_j |
1537 |
real( kind = dp) :: rho_i, rho_j |
1538 |
real( kind = dp) :: fshift_i, fshift_j |
1539 |
real( kind = dp) :: p_vdw, p_elect, p_hb, p_met |
1540 |
integer :: atid_i, atid_j, id1, id2, idx |
1541 |
integer :: k |
1542 |
|
1543 |
integer :: iHash |
1544 |
|
1545 |
r = sqrt(rijsq) |
1546 |
|
1547 |
vpair = 0.0_dp |
1548 |
fpair(1:3) = 0.0_dp |
1549 |
|
1550 |
p_vdw = 0.0 |
1551 |
p_elect = 0.0 |
1552 |
p_hb = 0.0 |
1553 |
p_met = 0.0 |
1554 |
|
1555 |
f1(1:3) = 0.0 |
1556 |
t1(1:3) = 0.0 |
1557 |
t2(1:3) = 0.0 |
1558 |
|
1559 |
#ifdef IS_MPI |
1560 |
atid_i = atid_row(i) |
1561 |
atid_j = atid_col(j) |
1562 |
|
1563 |
do idx = 1, 9 |
1564 |
A1(idx) = A_Row(idx, i) |
1565 |
A2(idx) = A_Col(idx, j) |
1566 |
eF1(idx) = eFrame_Row(idx, i) |
1567 |
eF2(idx) = eFrame_Col(idx, j) |
1568 |
enddo |
1569 |
|
1570 |
#else |
1571 |
atid_i = atid(i) |
1572 |
atid_j = atid(j) |
1573 |
do idx = 1, 9 |
1574 |
A1(idx) = A(idx, i) |
1575 |
A2(idx) = A(idx, j) |
1576 |
eF1(idx) = eFrame(idx, i) |
1577 |
eF2(idx) = eFrame(idx, j) |
1578 |
enddo |
1579 |
|
1580 |
#endif |
1581 |
|
1582 |
|
1583 |
iHash = InteractionHash(atid_i, atid_j) |
1584 |
|
1585 |
!! For the metallic potentials, we need to pass dF[rho]/drho since |
1586 |
!! the pair calculation routines no longer are aware of parallel. |
1587 |
|
1588 |
if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then |
1589 |
#ifdef IS_MPI |
1590 |
dfrhodrho_i = dfrhodrho_row(i) |
1591 |
dfrhodrho_j = dfrhodrho_col(j) |
1592 |
rho_i = rho_row(i) |
1593 |
rho_j = rho_col(j) |
1594 |
#else |
1595 |
dfrhodrho_i = dfrhodrho(i) |
1596 |
dfrhodrho_j = dfrhodrho(j) |
1597 |
rho_i = rho(i) |
1598 |
rho_j = rho(j) |
1599 |
#endif |
1600 |
end if |
1601 |
|
1602 |
vdwMult = vdwScale(topoDist) |
1603 |
electroMult = electrostaticScale(topoDist) |
1604 |
|
1605 |
if ( iand(iHash, LJ_PAIR).ne.0 ) then |
1606 |
call do_lj_pair(atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, & |
1607 |
p_vdw, f1) |
1608 |
endif |
1609 |
|
1610 |
if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then |
1611 |
call doElectrostaticPair(atid_i, atid_j, d, r, rijsq, rcut, sw, electroMult, & |
1612 |
vpair, p_elect, eF1, eF2, f1, t1, t2) |
1613 |
endif |
1614 |
|
1615 |
if ( iand(iHash, STICKY_PAIR).ne.0 ) then |
1616 |
call do_sticky_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, & |
1617 |
p_hb, A1, A2, f1, t1, t2) |
1618 |
endif |
1619 |
|
1620 |
if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then |
1621 |
call do_sticky_power_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, & |
1622 |
p_hb, A1, A2, f1, t1, t2) |
1623 |
endif |
1624 |
|
1625 |
if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then |
1626 |
call do_gb_pair(atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, & |
1627 |
p_vdw, A1, A2, f1, t1, t2) |
1628 |
endif |
1629 |
|
1630 |
if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then |
1631 |
call do_gb_pair(atid_i, atid_j, d, r, rijsq, sw, vdwMult, vpair, & |
1632 |
p_vdw, A1, A2, f1, t1, t2) |
1633 |
endif |
1634 |
|
1635 |
if ( iand(iHash, SHAPE_PAIR).ne.0 ) then |
1636 |
call do_shape_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, & |
1637 |
p_vdw, A1, A2, f1, t1, t2) |
1638 |
endif |
1639 |
|
1640 |
if ( iand(iHash, SHAPE_LJ).ne.0 ) then |
1641 |
call do_shape_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, & |
1642 |
p_vdw, A1, A2, f1, t1, t2) |
1643 |
endif |
1644 |
|
1645 |
if ( iand(iHash, EAM_PAIR).ne.0 ) then |
1646 |
call do_eam_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, & |
1647 |
p_met, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, fshift_i,fshift_j) |
1648 |
endif |
1649 |
|
1650 |
if ( iand(iHash, SC_PAIR).ne.0 ) then |
1651 |
call do_SC_pair(atid_i, atid_j, d, r, rijsq, sw, vpair, & |
1652 |
p_met, f1, rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, fshift_i, fshift_j) |
1653 |
endif |
1654 |
|
1655 |
if ( iand(iHash, MNM_PAIR).ne.0 ) then |
1656 |
call do_mnm_pair(atid_i, atid_j, d, r, rijsq, rcut, sw, vdwMult, vpair, & |
1657 |
p_vdw, A1, A2, f1, t1, t2) |
1658 |
endif |
1659 |
|
1660 |
|
1661 |
#ifdef IS_MPI |
1662 |
id1 = AtomRowToGlobal(i) |
1663 |
id2 = AtomColToGlobal(j) |
1664 |
|
1665 |
pot_row(VDW_POT,i) = pot_row(VDW_POT,i) + 0.5*p_vdw |
1666 |
pot_col(VDW_POT,j) = pot_col(VDW_POT,j) + 0.5*p_vdw |
1667 |
pot_row(ELECTROSTATIC_POT,i) = pot_row(ELECTROSTATIC_POT,i) + 0.5*p_elect |
1668 |
pot_col(ELECTROSTATIC_POT,j) = pot_col(ELECTROSTATIC_POT,j) + 0.5*p_elect |
1669 |
pot_row(HB_POT,i) = pot_row(HB_POT,i) + 0.5*p_hb |
1670 |
pot_col(HB_POT,j) = pot_col(HB_POT,j) + 0.5*p_hb |
1671 |
pot_Row(METALLIC_POT,i) = pot_Row(METALLIC_POT,i) + 0.5*p_met |
1672 |
pot_Col(METALLIC_POT,j) = pot_Col(METALLIC_POT,j) + 0.5*p_met |
1673 |
|
1674 |
do idx = 1, 3 |
1675 |
f_Row(idx,i) = f_Row(idx,i) + f1(idx) |
1676 |
f_Col(idx,j) = f_Col(idx,j) - f1(idx) |
1677 |
|
1678 |
t_Row(idx,i) = t_Row(idx,i) + t1(idx) |
1679 |
t_Col(idx,j) = t_Col(idx,j) + t2(idx) |
1680 |
enddo |
1681 |
! particle_pot is the difference between the full potential |
1682 |
! and the full potential without the presence of a particular |
1683 |
! particle (atom1). |
1684 |
! |
1685 |
! This reduces the density at other particle locations, so |
1686 |
! we need to recompute the density at atom2 assuming atom1 |
1687 |
! didn't contribute. This then requires recomputing the |
1688 |
! density functional for atom2 as well. |
1689 |
! |
1690 |
! Most of the particle_pot heavy lifting comes from the |
1691 |
! pair interaction, and will be handled by vpair. Parallel version. |
1692 |
|
1693 |
if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then |
1694 |
ppot_row(i) = ppot_row(i) - frho_row(j) + fshift_j |
1695 |
ppot_col(j) = ppot_col(j) - frho_col(i) + fshift_i |
1696 |
end if |
1697 |
|
1698 |
#else |
1699 |
id1 = i |
1700 |
id2 = j |
1701 |
|
1702 |
pot(VDW_POT) = pot(VDW_POT) + p_vdw |
1703 |
pot(ELECTROSTATIC_POT) = pot(ELECTROSTATIC_POT) + p_elect |
1704 |
pot(HB_POT) = pot(HB_POT) + p_hb |
1705 |
pot(METALLIC_POT) = pot(METALLIC_POT) + p_met |
1706 |
|
1707 |
do idx = 1, 3 |
1708 |
f(idx,i) = f(idx,i) + f1(idx) |
1709 |
f(idx,j) = f(idx,j) - f1(idx) |
1710 |
|
1711 |
t(idx,i) = t(idx,i) + t1(idx) |
1712 |
t(idx,j) = t(idx,j) + t2(idx) |
1713 |
enddo |
1714 |
! particle_pot is the difference between the full potential |
1715 |
! and the full potential without the presence of a particular |
1716 |
! particle (atom1). |
1717 |
! |
1718 |
! This reduces the density at other particle locations, so |
1719 |
! we need to recompute the density at atom2 assuming atom1 |
1720 |
! didn't contribute. This then requires recomputing the |
1721 |
! density functional for atom2 as well. |
1722 |
! |
1723 |
! Most of the particle_pot heavy lifting comes from the |
1724 |
! pair interaction, and will be handled by vpair. NonParallel version. |
1725 |
|
1726 |
if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then |
1727 |
particle_pot(i) = particle_pot(i) - frho(j) + fshift_j |
1728 |
particle_pot(j) = particle_pot(j) - frho(i) + fshift_i |
1729 |
end if |
1730 |
|
1731 |
|
1732 |
#endif |
1733 |
|
1734 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
1735 |
|
1736 |
fpair(1) = fpair(1) + f1(1) |
1737 |
fpair(2) = fpair(2) + f1(2) |
1738 |
fpair(3) = fpair(3) + f1(3) |
1739 |
|
1740 |
endif |
1741 |
|
1742 |
|
1743 |
!!$ |
1744 |
!!$ particle_pot(i) = particle_pot(i) + vpair*sw |
1745 |
!!$ particle_pot(j) = particle_pot(j) + vpair*sw |
1746 |
|
1747 |
end subroutine do_pair |
1748 |
|
1749 |
subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, & |
1750 |
eFrame, A, f, t, pot) |
1751 |
|
1752 |
real( kind = dp ) :: sw |
1753 |
real( kind = dp ), dimension(LR_POT_TYPES) :: pot |
1754 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
1755 |
real (kind=dp), dimension(9,nLocal) :: A |
1756 |
real (kind=dp), dimension(3,nLocal) :: f |
1757 |
real (kind=dp), dimension(3,nLocal) :: t |
1758 |
|
1759 |
integer, intent(in) :: i, j |
1760 |
real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut |
1761 |
real ( kind = dp ) :: r, rc |
1762 |
real ( kind = dp ), intent(inout) :: d(3), dc(3) |
1763 |
real ( kind = dp ) :: rho_i_at_j, rho_j_at_i |
1764 |
integer :: atid_i, atid_j, iHash |
1765 |
|
1766 |
r = sqrt(rijsq) |
1767 |
|
1768 |
#ifdef IS_MPI |
1769 |
atid_i = atid_row(i) |
1770 |
atid_j = atid_col(j) |
1771 |
#else |
1772 |
atid_i = atid(i) |
1773 |
atid_j = atid(j) |
1774 |
#endif |
1775 |
rho_i_at_j = 0.0_dp |
1776 |
rho_j_at_i = 0.0_dp |
1777 |
|
1778 |
iHash = InteractionHash(atid_i, atid_j) |
1779 |
|
1780 |
if ( iand(iHash, EAM_PAIR).ne.0 ) then |
1781 |
call calc_EAM_prepair_rho(atid_i, atid_j, d, r, rijsq, rho_i_at_j, rho_j_at_i) |
1782 |
endif |
1783 |
|
1784 |
if ( iand(iHash, SC_PAIR).ne.0 ) then |
1785 |
call calc_SC_prepair_rho(atid_i, atid_j, d, r, rijsq, rho_i_at_j, rho_j_at_i) |
1786 |
endif |
1787 |
|
1788 |
if ( iand(iHash, EAM_PAIR).ne.0 .or. iand(iHash, SC_PAIR).ne.0 ) then |
1789 |
#ifdef IS_MPI |
1790 |
rho_col(j) = rho_col(j) + rho_i_at_j |
1791 |
rho_row(i) = rho_row(i) + rho_j_at_i |
1792 |
#else |
1793 |
rho(j) = rho(j) + rho_i_at_j |
1794 |
rho(i) = rho(i) + rho_j_at_i |
1795 |
#endif |
1796 |
endif |
1797 |
|
1798 |
end subroutine do_prepair |
1799 |
|
1800 |
|
1801 |
subroutine do_preforce(nlocal, pot, particle_pot) |
1802 |
integer :: nlocal |
1803 |
real( kind = dp ),dimension(LR_POT_TYPES) :: pot |
1804 |
real( kind = dp ),dimension(nlocal) :: particle_pot |
1805 |
integer :: sc_err = 0 |
1806 |
|
1807 |
#ifdef IS_MPI |
1808 |
if ((FF_uses_EAM .and. SIM_uses_EAM) .or. (FF_uses_SC .and. SIM_uses_SC)) then |
1809 |
call scatter(rho_row,rho,plan_atom_row,sc_err) |
1810 |
if (sc_err /= 0 ) then |
1811 |
call handleError("do_preforce()", "Error scattering rho_row into rho") |
1812 |
endif |
1813 |
call scatter(rho_col,rho_tmp,plan_atom_col,sc_err) |
1814 |
if (sc_err /= 0 ) then |
1815 |
call handleError("do_preforce()", "Error scattering rho_col into rho") |
1816 |
endif |
1817 |
rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal) |
1818 |
end if |
1819 |
#endif |
1820 |
|
1821 |
|
1822 |
|
1823 |
if (FF_uses_EAM .and. SIM_uses_EAM) then |
1824 |
call calc_EAM_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot) |
1825 |
endif |
1826 |
if (FF_uses_SC .and. SIM_uses_SC) then |
1827 |
call calc_SC_preforce_Frho(nlocal, pot(METALLIC_POT), particle_pot) |
1828 |
endif |
1829 |
|
1830 |
#ifdef IS_MPI |
1831 |
if ((FF_uses_EAM .and. SIM_uses_EAM) .or. (FF_uses_SC .and. SIM_uses_SC)) then |
1832 |
!! communicate f(rho) and derivatives back into row and column arrays |
1833 |
call gather(frho,frho_row,plan_atom_row, sc_err) |
1834 |
if (sc_err /= 0) then |
1835 |
call handleError("do_preforce()","MPI gather frho_row failure") |
1836 |
endif |
1837 |
call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err) |
1838 |
if (sc_err /= 0) then |
1839 |
call handleError("do_preforce()","MPI gather dfrhodrho_row failure") |
1840 |
endif |
1841 |
call gather(frho,frho_col,plan_atom_col, sc_err) |
1842 |
if (sc_err /= 0) then |
1843 |
call handleError("do_preforce()","MPI gather frho_col failure") |
1844 |
endif |
1845 |
call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err) |
1846 |
if (sc_err /= 0) then |
1847 |
call handleError("do_preforce()","MPI gather dfrhodrho_col failure") |
1848 |
endif |
1849 |
end if |
1850 |
#endif |
1851 |
|
1852 |
end subroutine do_preforce |
1853 |
|
1854 |
|
1855 |
subroutine get_interatomic_vector(q_i, q_j, d, r_sq) |
1856 |
|
1857 |
real (kind = dp), dimension(3) :: q_i |
1858 |
real (kind = dp), dimension(3) :: q_j |
1859 |
real ( kind = dp ), intent(out) :: r_sq |
1860 |
real( kind = dp ) :: d(3), scaled(3) |
1861 |
integer i |
1862 |
|
1863 |
d(1) = q_j(1) - q_i(1) |
1864 |
d(2) = q_j(2) - q_i(2) |
1865 |
d(3) = q_j(3) - q_i(3) |
1866 |
|
1867 |
! Wrap back into periodic box if necessary |
1868 |
if ( SIM_uses_PBC ) then |
1869 |
|
1870 |
if( .not.boxIsOrthorhombic ) then |
1871 |
! calc the scaled coordinates. |
1872 |
! scaled = matmul(HmatInv, d) |
1873 |
|
1874 |
scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3) |
1875 |
scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3) |
1876 |
scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3) |
1877 |
|
1878 |
! wrap the scaled coordinates |
1879 |
|
1880 |
scaled(1) = scaled(1) - anint(scaled(1), kind=dp) |
1881 |
scaled(2) = scaled(2) - anint(scaled(2), kind=dp) |
1882 |
scaled(3) = scaled(3) - anint(scaled(3), kind=dp) |
1883 |
|
1884 |
! calc the wrapped real coordinates from the wrapped scaled |
1885 |
! coordinates |
1886 |
! d = matmul(Hmat,scaled) |
1887 |
d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3) |
1888 |
d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3) |
1889 |
d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3) |
1890 |
|
1891 |
else |
1892 |
! calc the scaled coordinates. |
1893 |
|
1894 |
scaled(1) = d(1) * HmatInv(1,1) |
1895 |
scaled(2) = d(2) * HmatInv(2,2) |
1896 |
scaled(3) = d(3) * HmatInv(3,3) |
1897 |
|
1898 |
! wrap the scaled coordinates |
1899 |
|
1900 |
scaled(1) = scaled(1) - anint(scaled(1), kind=dp) |
1901 |
scaled(2) = scaled(2) - anint(scaled(2), kind=dp) |
1902 |
scaled(3) = scaled(3) - anint(scaled(3), kind=dp) |
1903 |
|
1904 |
! calc the wrapped real coordinates from the wrapped scaled |
1905 |
! coordinates |
1906 |
|
1907 |
d(1) = scaled(1)*Hmat(1,1) |
1908 |
d(2) = scaled(2)*Hmat(2,2) |
1909 |
d(3) = scaled(3)*Hmat(3,3) |
1910 |
|
1911 |
endif |
1912 |
|
1913 |
endif |
1914 |
|
1915 |
r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3) |
1916 |
|
1917 |
end subroutine get_interatomic_vector |
1918 |
|
1919 |
subroutine zero_work_arrays() |
1920 |
|
1921 |
#ifdef IS_MPI |
1922 |
|
1923 |
q_Row = 0.0_dp |
1924 |
q_Col = 0.0_dp |
1925 |
|
1926 |
q_group_Row = 0.0_dp |
1927 |
q_group_Col = 0.0_dp |
1928 |
|
1929 |
eFrame_Row = 0.0_dp |
1930 |
eFrame_Col = 0.0_dp |
1931 |
|
1932 |
A_Row = 0.0_dp |
1933 |
A_Col = 0.0_dp |
1934 |
|
1935 |
f_Row = 0.0_dp |
1936 |
f_Col = 0.0_dp |
1937 |
f_Temp = 0.0_dp |
1938 |
|
1939 |
t_Row = 0.0_dp |
1940 |
t_Col = 0.0_dp |
1941 |
t_Temp = 0.0_dp |
1942 |
|
1943 |
pot_Row = 0.0_dp |
1944 |
pot_Col = 0.0_dp |
1945 |
pot_Temp = 0.0_dp |
1946 |
ppot_Temp = 0.0_dp |
1947 |
|
1948 |
frho_row = 0.0_dp |
1949 |
frho_col = 0.0_dp |
1950 |
rho_row = 0.0_dp |
1951 |
rho_col = 0.0_dp |
1952 |
rho_tmp = 0.0_dp |
1953 |
dfrhodrho_row = 0.0_dp |
1954 |
dfrhodrho_col = 0.0_dp |
1955 |
|
1956 |
#endif |
1957 |
rho = 0.0_dp |
1958 |
frho = 0.0_dp |
1959 |
dfrhodrho = 0.0_dp |
1960 |
|
1961 |
end subroutine zero_work_arrays |
1962 |
|
1963 |
function skipThisPair(atom1, atom2) result(skip_it) |
1964 |
integer, intent(in) :: atom1 |
1965 |
integer, intent(in), optional :: atom2 |
1966 |
logical :: skip_it |
1967 |
integer :: unique_id_1, unique_id_2 |
1968 |
integer :: me_i,me_j |
1969 |
integer :: i |
1970 |
|
1971 |
skip_it = .false. |
1972 |
|
1973 |
!! there are a number of reasons to skip a pair or a particle |
1974 |
!! mostly we do this to exclude atoms who are involved in short |
1975 |
!! range interactions (bonds, bends, torsions), but we also need |
1976 |
!! to exclude some overcounted interactions that result from |
1977 |
!! the parallel decomposition |
1978 |
|
1979 |
#ifdef IS_MPI |
1980 |
!! in MPI, we have to look up the unique IDs for each atom |
1981 |
unique_id_1 = AtomRowToGlobal(atom1) |
1982 |
unique_id_2 = AtomColToGlobal(atom2) |
1983 |
!! this situation should only arise in MPI simulations |
1984 |
if (unique_id_1 == unique_id_2) then |
1985 |
skip_it = .true. |
1986 |
return |
1987 |
end if |
1988 |
|
1989 |
!! this prevents us from doing the pair on multiple processors |
1990 |
if (unique_id_1 < unique_id_2) then |
1991 |
if (mod(unique_id_1 + unique_id_2,2) == 0) then |
1992 |
skip_it = .true. |
1993 |
return |
1994 |
endif |
1995 |
else |
1996 |
if (mod(unique_id_1 + unique_id_2,2) == 1) then |
1997 |
skip_it = .true. |
1998 |
return |
1999 |
endif |
2000 |
endif |
2001 |
#else |
2002 |
!! in the normal loop, the atom numbers are unique |
2003 |
unique_id_1 = atom1 |
2004 |
unique_id_2 = atom2 |
2005 |
#endif |
2006 |
|
2007 |
#ifdef IS_MPI |
2008 |
do i = 1, nSkipsForRowAtom(atom1) |
2009 |
if (skipsForRowAtom(atom1, i) .eq. unique_id_2) then |
2010 |
skip_it = .true. |
2011 |
return |
2012 |
endif |
2013 |
end do |
2014 |
#else |
2015 |
do i = 1, nSkipsForLocalAtom(atom1) |
2016 |
if (skipsForLocalAtom(atom1, i) .eq. unique_id_2) then |
2017 |
skip_it = .true. |
2018 |
return |
2019 |
endif |
2020 |
end do |
2021 |
#endif |
2022 |
|
2023 |
return |
2024 |
end function skipThisPair |
2025 |
|
2026 |
function getTopoDistance(atom1, atom2) result(topoDist) |
2027 |
integer, intent(in) :: atom1 |
2028 |
integer, intent(in) :: atom2 |
2029 |
integer :: topoDist |
2030 |
integer :: unique_id_2 |
2031 |
integer :: i |
2032 |
|
2033 |
#ifdef IS_MPI |
2034 |
unique_id_2 = AtomColToGlobal(atom2) |
2035 |
#else |
2036 |
unique_id_2 = atom2 |
2037 |
#endif |
2038 |
|
2039 |
! zero is default for unconnected (i.e. normal) pair interactions |
2040 |
|
2041 |
topoDist = 0 |
2042 |
|
2043 |
do i = 1, nTopoPairsForAtom(atom1) |
2044 |
if (toposForAtom(atom1, i) .eq. unique_id_2) then |
2045 |
topoDist = topoDistance(atom1, i) |
2046 |
return |
2047 |
endif |
2048 |
end do |
2049 |
|
2050 |
return |
2051 |
end function getTopoDistance |
2052 |
|
2053 |
function FF_UsesDirectionalAtoms() result(doesit) |
2054 |
logical :: doesit |
2055 |
doesit = FF_uses_DirectionalAtoms |
2056 |
end function FF_UsesDirectionalAtoms |
2057 |
|
2058 |
function FF_RequiresPrepairCalc() result(doesit) |
2059 |
logical :: doesit |
2060 |
doesit = FF_uses_EAM .or. FF_uses_SC |
2061 |
end function FF_RequiresPrepairCalc |
2062 |
|
2063 |
#ifdef PROFILE |
2064 |
function getforcetime() result(totalforcetime) |
2065 |
real(kind=dp) :: totalforcetime |
2066 |
totalforcetime = forcetime |
2067 |
end function getforcetime |
2068 |
#endif |
2069 |
|
2070 |
!! This cleans componets of force arrays belonging only to fortran |
2071 |
|
2072 |
subroutine add_stress_tensor(dpair, fpair, tau) |
2073 |
|
2074 |
real( kind = dp ), dimension(3), intent(in) :: dpair, fpair |
2075 |
real( kind = dp ), dimension(9), intent(inout) :: tau |
2076 |
|
2077 |
! because the d vector is the rj - ri vector, and |
2078 |
! because fx, fy, fz are the force on atom i, we need a |
2079 |
! negative sign here: |
2080 |
|
2081 |
tau(1) = tau(1) - dpair(1) * fpair(1) |
2082 |
tau(2) = tau(2) - dpair(1) * fpair(2) |
2083 |
tau(3) = tau(3) - dpair(1) * fpair(3) |
2084 |
tau(4) = tau(4) - dpair(2) * fpair(1) |
2085 |
tau(5) = tau(5) - dpair(2) * fpair(2) |
2086 |
tau(6) = tau(6) - dpair(2) * fpair(3) |
2087 |
tau(7) = tau(7) - dpair(3) * fpair(1) |
2088 |
tau(8) = tau(8) - dpair(3) * fpair(2) |
2089 |
tau(9) = tau(9) - dpair(3) * fpair(3) |
2090 |
|
2091 |
end subroutine add_stress_tensor |
2092 |
|
2093 |
end module doForces |