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