1 |
!! |
2 |
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
!! |
4 |
!! The University of Notre Dame grants you ("Licensee") a |
5 |
!! non-exclusive, royalty free, license to use, modify and |
6 |
!! redistribute this software in source and binary code form, provided |
7 |
!! that the following conditions are met: |
8 |
!! |
9 |
!! 1. Acknowledgement of the program authors must be made in any |
10 |
!! publication of scientific results based in part on use of the |
11 |
!! program. An acceptable form of acknowledgement is citation of |
12 |
!! the article in which the program was described (Matthew |
13 |
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
!! Parallel Simulation Engine for Molecular Dynamics," |
16 |
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
!! |
18 |
!! 2. Redistributions of source code must retain the above copyright |
19 |
!! notice, this list of conditions and the following disclaimer. |
20 |
!! |
21 |
!! 3. Redistributions in binary form must reproduce the above copyright |
22 |
!! notice, this list of conditions and the following disclaimer in the |
23 |
!! documentation and/or other materials provided with the |
24 |
!! distribution. |
25 |
!! |
26 |
!! This software is provided "AS IS," without a warranty of any |
27 |
!! kind. All express or implied conditions, representations and |
28 |
!! warranties, including any implied warranty of merchantability, |
29 |
!! fitness for a particular purpose or non-infringement, are hereby |
30 |
!! excluded. The University of Notre Dame and its licensors shall not |
31 |
!! be liable for any damages suffered by licensee as a result of |
32 |
!! using, modifying or distributing the software or its |
33 |
!! derivatives. In no event will the University of Notre Dame or its |
34 |
!! licensors be liable for any lost revenue, profit or data, or for |
35 |
!! direct, indirect, special, consequential, incidental or punitive |
36 |
!! damages, however caused and regardless of the theory of liability, |
37 |
!! arising out of the use of or inability to use software, even if the |
38 |
!! University of Notre Dame has been advised of the possibility of |
39 |
!! such damages. |
40 |
!! |
41 |
|
42 |
!! doForces.F90 |
43 |
!! module doForces |
44 |
!! Calculates Long Range forces. |
45 |
|
46 |
!! @author Charles F. Vardeman II |
47 |
!! @author Matthew Meineke |
48 |
!! @version $Id: doForces.F90,v 1.66 2005-11-02 20:35:34 chrisfen Exp $, $Date: 2005-11-02 20:35:34 $, $Name: not supported by cvs2svn $, $Revision: 1.66 $ |
49 |
|
50 |
|
51 |
module doForces |
52 |
use force_globals |
53 |
use simulation |
54 |
use definitions |
55 |
use atype_module |
56 |
use switcheroo |
57 |
use neighborLists |
58 |
use lj |
59 |
use sticky |
60 |
use electrostatic_module |
61 |
use gayberne |
62 |
use shapes |
63 |
use vector_class |
64 |
use eam |
65 |
use status |
66 |
#ifdef IS_MPI |
67 |
use mpiSimulation |
68 |
#endif |
69 |
|
70 |
implicit none |
71 |
PRIVATE |
72 |
|
73 |
#define __FORTRAN90 |
74 |
#include "UseTheForce/fSwitchingFunction.h" |
75 |
#include "UseTheForce/fCutoffPolicy.h" |
76 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
77 |
#include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h" |
78 |
|
79 |
|
80 |
INTEGER, PARAMETER:: PREPAIR_LOOP = 1 |
81 |
INTEGER, PARAMETER:: PAIR_LOOP = 2 |
82 |
|
83 |
logical, save :: haveNeighborList = .false. |
84 |
logical, save :: haveSIMvariables = .false. |
85 |
logical, save :: haveSaneForceField = .false. |
86 |
logical, save :: haveInteractionHash = .false. |
87 |
logical, save :: haveGtypeCutoffMap = .false. |
88 |
logical, save :: haveDefaultCutoffs = .false. |
89 |
logical, save :: haveRlist = .false. |
90 |
|
91 |
logical, save :: FF_uses_DirectionalAtoms |
92 |
logical, save :: FF_uses_Dipoles |
93 |
logical, save :: FF_uses_GayBerne |
94 |
logical, save :: FF_uses_EAM |
95 |
|
96 |
logical, save :: SIM_uses_DirectionalAtoms |
97 |
logical, save :: SIM_uses_EAM |
98 |
logical, save :: SIM_requires_postpair_calc |
99 |
logical, save :: SIM_requires_prepair_calc |
100 |
logical, save :: SIM_uses_PBC |
101 |
|
102 |
integer, save :: electrostaticSummationMethod |
103 |
|
104 |
public :: init_FF |
105 |
public :: setDefaultCutoffs |
106 |
public :: do_force_loop |
107 |
public :: createInteractionHash |
108 |
public :: createGtypeCutoffMap |
109 |
public :: getStickyCut |
110 |
public :: getStickyPowerCut |
111 |
public :: getGayBerneCut |
112 |
public :: getEAMCut |
113 |
public :: getShapeCut |
114 |
|
115 |
#ifdef PROFILE |
116 |
public :: getforcetime |
117 |
real, save :: forceTime = 0 |
118 |
real :: forceTimeInitial, forceTimeFinal |
119 |
integer :: nLoops |
120 |
#endif |
121 |
|
122 |
!! Variables for cutoff mapping and interaction mapping |
123 |
! Bit hash to determine pair-pair interactions. |
124 |
integer, dimension(:,:), allocatable :: InteractionHash |
125 |
real(kind=dp), dimension(:), allocatable :: atypeMaxCutoff |
126 |
real(kind=dp), dimension(:), allocatable, target :: groupMaxCutoffRow |
127 |
real(kind=dp), dimension(:), pointer :: groupMaxCutoffCol |
128 |
|
129 |
integer, dimension(:), allocatable, target :: groupToGtypeRow |
130 |
integer, dimension(:), pointer :: groupToGtypeCol => null() |
131 |
|
132 |
real(kind=dp), dimension(:), allocatable,target :: gtypeMaxCutoffRow |
133 |
real(kind=dp), dimension(:), pointer :: gtypeMaxCutoffCol |
134 |
type ::gtypeCutoffs |
135 |
real(kind=dp) :: rcut |
136 |
real(kind=dp) :: rcutsq |
137 |
real(kind=dp) :: rlistsq |
138 |
end type gtypeCutoffs |
139 |
type(gtypeCutoffs), dimension(:,:), allocatable :: gtypeCutoffMap |
140 |
|
141 |
integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY |
142 |
real(kind=dp),save :: defaultRcut, defaultRsw, defaultRlist |
143 |
real(kind=dp),save :: listSkin |
144 |
|
145 |
contains |
146 |
|
147 |
subroutine createInteractionHash(status) |
148 |
integer :: nAtypes |
149 |
integer, intent(out) :: status |
150 |
integer :: i |
151 |
integer :: j |
152 |
integer :: iHash |
153 |
!! Test Types |
154 |
logical :: i_is_LJ |
155 |
logical :: i_is_Elect |
156 |
logical :: i_is_Sticky |
157 |
logical :: i_is_StickyP |
158 |
logical :: i_is_GB |
159 |
logical :: i_is_EAM |
160 |
logical :: i_is_Shape |
161 |
logical :: j_is_LJ |
162 |
logical :: j_is_Elect |
163 |
logical :: j_is_Sticky |
164 |
logical :: j_is_StickyP |
165 |
logical :: j_is_GB |
166 |
logical :: j_is_EAM |
167 |
logical :: j_is_Shape |
168 |
real(kind=dp) :: myRcut |
169 |
|
170 |
status = 0 |
171 |
|
172 |
if (.not. associated(atypes)) then |
173 |
call handleError("atype", "atypes was not present before call of createInteractionHash!") |
174 |
status = -1 |
175 |
return |
176 |
endif |
177 |
|
178 |
nAtypes = getSize(atypes) |
179 |
|
180 |
if (nAtypes == 0) then |
181 |
status = -1 |
182 |
return |
183 |
end if |
184 |
|
185 |
if (.not. allocated(InteractionHash)) then |
186 |
allocate(InteractionHash(nAtypes,nAtypes)) |
187 |
else |
188 |
deallocate(InteractionHash) |
189 |
allocate(InteractionHash(nAtypes,nAtypes)) |
190 |
endif |
191 |
|
192 |
if (.not. allocated(atypeMaxCutoff)) then |
193 |
allocate(atypeMaxCutoff(nAtypes)) |
194 |
else |
195 |
deallocate(atypeMaxCutoff) |
196 |
allocate(atypeMaxCutoff(nAtypes)) |
197 |
endif |
198 |
|
199 |
do i = 1, nAtypes |
200 |
call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ) |
201 |
call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect) |
202 |
call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky) |
203 |
call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP) |
204 |
call getElementProperty(atypes, i, "is_GayBerne", i_is_GB) |
205 |
call getElementProperty(atypes, i, "is_EAM", i_is_EAM) |
206 |
call getElementProperty(atypes, i, "is_Shape", i_is_Shape) |
207 |
|
208 |
do j = i, nAtypes |
209 |
|
210 |
iHash = 0 |
211 |
myRcut = 0.0_dp |
212 |
|
213 |
call getElementProperty(atypes, j, "is_LennardJones", j_is_LJ) |
214 |
call getElementProperty(atypes, j, "is_Electrostatic", j_is_Elect) |
215 |
call getElementProperty(atypes, j, "is_Sticky", j_is_Sticky) |
216 |
call getElementProperty(atypes, j, "is_StickyPower", j_is_StickyP) |
217 |
call getElementProperty(atypes, j, "is_GayBerne", j_is_GB) |
218 |
call getElementProperty(atypes, j, "is_EAM", j_is_EAM) |
219 |
call getElementProperty(atypes, j, "is_Shape", j_is_Shape) |
220 |
|
221 |
if (i_is_LJ .and. j_is_LJ) then |
222 |
iHash = ior(iHash, LJ_PAIR) |
223 |
endif |
224 |
|
225 |
if (i_is_Elect .and. j_is_Elect) then |
226 |
iHash = ior(iHash, ELECTROSTATIC_PAIR) |
227 |
endif |
228 |
|
229 |
if (i_is_Sticky .and. j_is_Sticky) then |
230 |
iHash = ior(iHash, STICKY_PAIR) |
231 |
endif |
232 |
|
233 |
if (i_is_StickyP .and. j_is_StickyP) then |
234 |
iHash = ior(iHash, STICKYPOWER_PAIR) |
235 |
endif |
236 |
|
237 |
if (i_is_EAM .and. j_is_EAM) then |
238 |
iHash = ior(iHash, EAM_PAIR) |
239 |
endif |
240 |
|
241 |
if (i_is_GB .and. j_is_GB) iHash = ior(iHash, GAYBERNE_PAIR) |
242 |
if (i_is_GB .and. j_is_LJ) iHash = ior(iHash, GAYBERNE_LJ) |
243 |
if (i_is_LJ .and. j_is_GB) iHash = ior(iHash, GAYBERNE_LJ) |
244 |
|
245 |
if (i_is_Shape .and. j_is_Shape) iHash = ior(iHash, SHAPE_PAIR) |
246 |
if (i_is_Shape .and. j_is_LJ) iHash = ior(iHash, SHAPE_LJ) |
247 |
if (i_is_LJ .and. j_is_Shape) iHash = ior(iHash, SHAPE_LJ) |
248 |
|
249 |
|
250 |
InteractionHash(i,j) = iHash |
251 |
InteractionHash(j,i) = iHash |
252 |
|
253 |
end do |
254 |
|
255 |
end do |
256 |
|
257 |
haveInteractionHash = .true. |
258 |
end subroutine createInteractionHash |
259 |
|
260 |
subroutine createGtypeCutoffMap(stat) |
261 |
|
262 |
integer, intent(out), optional :: stat |
263 |
logical :: i_is_LJ |
264 |
logical :: i_is_Elect |
265 |
logical :: i_is_Sticky |
266 |
logical :: i_is_StickyP |
267 |
logical :: i_is_GB |
268 |
logical :: i_is_EAM |
269 |
logical :: i_is_Shape |
270 |
logical :: GtypeFound |
271 |
|
272 |
integer :: myStatus, nAtypes, i, j, istart, iend, jstart, jend |
273 |
integer :: n_in_i, me_i, ia, g, atom1, ja, n_in_j,me_j |
274 |
integer :: nGroupsInRow |
275 |
integer :: nGroupsInCol |
276 |
integer :: nGroupTypesRow,nGroupTypesCol |
277 |
real(kind=dp):: thisSigma, bigSigma, thisRcut, tradRcut, tol, skin |
278 |
real(kind=dp) :: biggestAtypeCutoff |
279 |
|
280 |
stat = 0 |
281 |
if (.not. haveInteractionHash) then |
282 |
call createInteractionHash(myStatus) |
283 |
if (myStatus .ne. 0) then |
284 |
write(default_error, *) 'createInteractionHash failed in doForces!' |
285 |
stat = -1 |
286 |
return |
287 |
endif |
288 |
endif |
289 |
#ifdef IS_MPI |
290 |
nGroupsInRow = getNgroupsInRow(plan_group_row) |
291 |
nGroupsInCol = getNgroupsInCol(plan_group_col) |
292 |
#endif |
293 |
nAtypes = getSize(atypes) |
294 |
! Set all of the initial cutoffs to zero. |
295 |
atypeMaxCutoff = 0.0_dp |
296 |
do i = 1, nAtypes |
297 |
if (SimHasAtype(i)) then |
298 |
call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ) |
299 |
call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect) |
300 |
call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky) |
301 |
call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP) |
302 |
call getElementProperty(atypes, i, "is_GayBerne", i_is_GB) |
303 |
call getElementProperty(atypes, i, "is_EAM", i_is_EAM) |
304 |
call getElementProperty(atypes, i, "is_Shape", i_is_Shape) |
305 |
|
306 |
|
307 |
if (haveDefaultCutoffs) then |
308 |
atypeMaxCutoff(i) = defaultRcut |
309 |
else |
310 |
if (i_is_LJ) then |
311 |
thisRcut = getSigma(i) * 2.5_dp |
312 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
313 |
endif |
314 |
if (i_is_Elect) then |
315 |
thisRcut = defaultRcut |
316 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
317 |
endif |
318 |
if (i_is_Sticky) then |
319 |
thisRcut = getStickyCut(i) |
320 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
321 |
endif |
322 |
if (i_is_StickyP) then |
323 |
thisRcut = getStickyPowerCut(i) |
324 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
325 |
endif |
326 |
if (i_is_GB) then |
327 |
thisRcut = getGayBerneCut(i) |
328 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
329 |
endif |
330 |
if (i_is_EAM) then |
331 |
thisRcut = getEAMCut(i) |
332 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
333 |
endif |
334 |
if (i_is_Shape) then |
335 |
thisRcut = getShapeCut(i) |
336 |
if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut |
337 |
endif |
338 |
endif |
339 |
|
340 |
|
341 |
if (atypeMaxCutoff(i).gt.biggestAtypeCutoff) then |
342 |
biggestAtypeCutoff = atypeMaxCutoff(i) |
343 |
endif |
344 |
|
345 |
endif |
346 |
enddo |
347 |
|
348 |
|
349 |
|
350 |
istart = 1 |
351 |
jstart = 1 |
352 |
#ifdef IS_MPI |
353 |
iend = nGroupsInRow |
354 |
jend = nGroupsInCol |
355 |
#else |
356 |
iend = nGroups |
357 |
jend = nGroups |
358 |
#endif |
359 |
|
360 |
!! allocate the groupToGtype and gtypeMaxCutoff here. |
361 |
if(.not.allocated(groupToGtypeRow)) then |
362 |
! allocate(groupToGtype(iend)) |
363 |
allocate(groupToGtypeRow(iend)) |
364 |
else |
365 |
deallocate(groupToGtypeRow) |
366 |
allocate(groupToGtypeRow(iend)) |
367 |
endif |
368 |
if(.not.allocated(groupMaxCutoffRow)) then |
369 |
allocate(groupMaxCutoffRow(iend)) |
370 |
else |
371 |
deallocate(groupMaxCutoffRow) |
372 |
allocate(groupMaxCutoffRow(iend)) |
373 |
end if |
374 |
|
375 |
if(.not.allocated(gtypeMaxCutoffRow)) then |
376 |
allocate(gtypeMaxCutoffRow(iend)) |
377 |
else |
378 |
deallocate(gtypeMaxCutoffRow) |
379 |
allocate(gtypeMaxCutoffRow(iend)) |
380 |
endif |
381 |
|
382 |
|
383 |
#ifdef IS_MPI |
384 |
! We only allocate new storage if we are in MPI because Ncol /= Nrow |
385 |
if(.not.associated(groupToGtypeCol)) then |
386 |
allocate(groupToGtypeCol(jend)) |
387 |
else |
388 |
deallocate(groupToGtypeCol) |
389 |
allocate(groupToGtypeCol(jend)) |
390 |
end if |
391 |
|
392 |
if(.not.associated(groupToGtypeCol)) then |
393 |
allocate(groupToGtypeCol(jend)) |
394 |
else |
395 |
deallocate(groupToGtypeCol) |
396 |
allocate(groupToGtypeCol(jend)) |
397 |
end if |
398 |
if(.not.associated(gtypeMaxCutoffCol)) then |
399 |
allocate(gtypeMaxCutoffCol(jend)) |
400 |
else |
401 |
deallocate(gtypeMaxCutoffCol) |
402 |
allocate(gtypeMaxCutoffCol(jend)) |
403 |
end if |
404 |
|
405 |
groupMaxCutoffCol = 0.0_dp |
406 |
gtypeMaxCutoffCol = 0.0_dp |
407 |
|
408 |
#endif |
409 |
groupMaxCutoffRow = 0.0_dp |
410 |
gtypeMaxCutoffRow = 0.0_dp |
411 |
|
412 |
|
413 |
!! first we do a single loop over the cutoff groups to find the |
414 |
!! largest cutoff for any atypes present in this group. We also |
415 |
!! create gtypes at this point. |
416 |
|
417 |
tol = 1.0d-6 |
418 |
nGroupTypesRow = 0 |
419 |
|
420 |
do i = istart, iend |
421 |
n_in_i = groupStartRow(i+1) - groupStartRow(i) |
422 |
groupMaxCutoffRow(i) = 0.0_dp |
423 |
do ia = groupStartRow(i), groupStartRow(i+1)-1 |
424 |
atom1 = groupListRow(ia) |
425 |
#ifdef IS_MPI |
426 |
me_i = atid_row(atom1) |
427 |
#else |
428 |
me_i = atid(atom1) |
429 |
#endif |
430 |
if (atypeMaxCutoff(me_i).gt.groupMaxCutoffRow(i)) then |
431 |
groupMaxCutoffRow(i)=atypeMaxCutoff(me_i) |
432 |
endif |
433 |
enddo |
434 |
|
435 |
if (nGroupTypesRow.eq.0) then |
436 |
nGroupTypesRow = nGroupTypesRow + 1 |
437 |
gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i) |
438 |
groupToGtypeRow(i) = nGroupTypesRow |
439 |
else |
440 |
GtypeFound = .false. |
441 |
do g = 1, nGroupTypesRow |
442 |
if ( abs(groupMaxCutoffRow(i) - gtypeMaxCutoffRow(g)).lt.tol) then |
443 |
groupToGtypeRow(i) = g |
444 |
GtypeFound = .true. |
445 |
endif |
446 |
enddo |
447 |
if (.not.GtypeFound) then |
448 |
nGroupTypesRow = nGroupTypesRow + 1 |
449 |
gtypeMaxCutoffRow(nGroupTypesRow) = groupMaxCutoffRow(i) |
450 |
groupToGtypeRow(i) = nGroupTypesRow |
451 |
endif |
452 |
endif |
453 |
enddo |
454 |
|
455 |
#ifdef IS_MPI |
456 |
do j = jstart, jend |
457 |
n_in_j = groupStartCol(j+1) - groupStartCol(j) |
458 |
groupMaxCutoffCol(j) = 0.0_dp |
459 |
do ja = groupStartCol(j), groupStartCol(j+1)-1 |
460 |
atom1 = groupListCol(ja) |
461 |
|
462 |
me_j = atid_col(atom1) |
463 |
|
464 |
if (atypeMaxCutoff(me_j).gt.groupMaxCutoffCol(j)) then |
465 |
groupMaxCutoffCol(j)=atypeMaxCutoff(me_j) |
466 |
endif |
467 |
enddo |
468 |
|
469 |
if (nGroupTypesCol.eq.0) then |
470 |
nGroupTypesCol = nGroupTypesCol + 1 |
471 |
gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j) |
472 |
groupToGtypeCol(j) = nGroupTypesCol |
473 |
else |
474 |
GtypeFound = .false. |
475 |
do g = 1, nGroupTypesCol |
476 |
if ( abs(groupMaxCutoffCol(j) - gtypeMaxCutoffCol(g)).lt.tol) then |
477 |
groupToGtypeCol(j) = g |
478 |
GtypeFound = .true. |
479 |
endif |
480 |
enddo |
481 |
if (.not.GtypeFound) then |
482 |
nGroupTypesCol = nGroupTypesCol + 1 |
483 |
gtypeMaxCutoffCol(nGroupTypesCol) = groupMaxCutoffCol(j) |
484 |
groupToGtypeCol(j) = nGroupTypesCol |
485 |
endif |
486 |
endif |
487 |
enddo |
488 |
|
489 |
#else |
490 |
! Set pointers to information we just found |
491 |
nGroupTypesCol = nGroupTypesRow |
492 |
groupToGtypeCol => groupToGtypeRow |
493 |
gtypeMaxCutoffCol => gtypeMaxCutoffRow |
494 |
groupMaxCutoffCol => groupMaxCutoffRow |
495 |
#endif |
496 |
|
497 |
|
498 |
|
499 |
|
500 |
|
501 |
!! allocate the gtypeCutoffMap here. |
502 |
allocate(gtypeCutoffMap(nGroupTypesRow,nGroupTypesCol)) |
503 |
!! then we do a double loop over all the group TYPES to find the cutoff |
504 |
!! map between groups of two types |
505 |
tradRcut = max(maxval(gtypeMaxCutoffRow),maxval(gtypeMaxCutoffCol)) |
506 |
|
507 |
do i = 1, nGroupTypesRow |
508 |
do j = 1, nGroupTypesCol |
509 |
|
510 |
select case(cutoffPolicy) |
511 |
case(TRADITIONAL_CUTOFF_POLICY) |
512 |
thisRcut = tradRcut |
513 |
case(MIX_CUTOFF_POLICY) |
514 |
thisRcut = 0.5_dp * (gtypeMaxCutoffRow(i) + gtypeMaxCutoffCol(j)) |
515 |
case(MAX_CUTOFF_POLICY) |
516 |
thisRcut = max(gtypeMaxCutoffRow(i), gtypeMaxCutoffCol(j)) |
517 |
case default |
518 |
call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy") |
519 |
return |
520 |
end select |
521 |
gtypeCutoffMap(i,j)%rcut = thisRcut |
522 |
gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut |
523 |
skin = defaultRlist - defaultRcut |
524 |
listSkin = skin ! set neighbor list skin thickness |
525 |
gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skin)**2 |
526 |
|
527 |
! sanity check |
528 |
|
529 |
if (haveDefaultCutoffs) then |
530 |
if (abs(gtypeCutoffMap(i,j)%rcut - defaultRcut).gt.0.0001) then |
531 |
call handleError("createGtypeCutoffMap", "user-specified rCut does not match computed group Cutoff") |
532 |
endif |
533 |
endif |
534 |
enddo |
535 |
enddo |
536 |
if(allocated(gtypeMaxCutoffRow)) deallocate(gtypeMaxCutoffRow) |
537 |
if(allocated(groupMaxCutoffRow)) deallocate(groupMaxCutoffRow) |
538 |
if(allocated(atypeMaxCutoff)) deallocate(atypeMaxCutoff) |
539 |
#ifdef IS_MPI |
540 |
if(associated(groupMaxCutoffCol)) deallocate(groupMaxCutoffCol) |
541 |
if(associated(gtypeMaxCutoffCol)) deallocate(gtypeMaxCutoffCol) |
542 |
#endif |
543 |
groupMaxCutoffCol => null() |
544 |
gtypeMaxCutoffCol => null() |
545 |
|
546 |
haveGtypeCutoffMap = .true. |
547 |
end subroutine createGtypeCutoffMap |
548 |
|
549 |
subroutine setDefaultCutoffs(defRcut, defRsw, defRlist, cutPolicy) |
550 |
real(kind=dp),intent(in) :: defRcut, defRsw, defRlist |
551 |
integer, intent(in) :: cutPolicy |
552 |
|
553 |
defaultRcut = defRcut |
554 |
defaultRsw = defRsw |
555 |
defaultRlist = defRlist |
556 |
cutoffPolicy = cutPolicy |
557 |
|
558 |
haveDefaultCutoffs = .true. |
559 |
end subroutine setDefaultCutoffs |
560 |
|
561 |
subroutine setCutoffPolicy(cutPolicy) |
562 |
|
563 |
integer, intent(in) :: cutPolicy |
564 |
cutoffPolicy = cutPolicy |
565 |
call createGtypeCutoffMap() |
566 |
end subroutine setCutoffPolicy |
567 |
|
568 |
|
569 |
subroutine setSimVariables() |
570 |
SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms() |
571 |
SIM_uses_EAM = SimUsesEAM() |
572 |
SIM_requires_postpair_calc = SimRequiresPostpairCalc() |
573 |
SIM_requires_prepair_calc = SimRequiresPrepairCalc() |
574 |
SIM_uses_PBC = SimUsesPBC() |
575 |
|
576 |
haveSIMvariables = .true. |
577 |
|
578 |
return |
579 |
end subroutine setSimVariables |
580 |
|
581 |
subroutine doReadyCheck(error) |
582 |
integer, intent(out) :: error |
583 |
|
584 |
integer :: myStatus |
585 |
|
586 |
error = 0 |
587 |
|
588 |
if (.not. haveInteractionHash) then |
589 |
myStatus = 0 |
590 |
call createInteractionHash(myStatus) |
591 |
if (myStatus .ne. 0) then |
592 |
write(default_error, *) 'createInteractionHash failed in doForces!' |
593 |
error = -1 |
594 |
return |
595 |
endif |
596 |
endif |
597 |
|
598 |
if (.not. haveGtypeCutoffMap) then |
599 |
myStatus = 0 |
600 |
call createGtypeCutoffMap(myStatus) |
601 |
if (myStatus .ne. 0) then |
602 |
write(default_error, *) 'createGtypeCutoffMap failed in doForces!' |
603 |
error = -1 |
604 |
return |
605 |
endif |
606 |
endif |
607 |
|
608 |
if (.not. haveSIMvariables) then |
609 |
call setSimVariables() |
610 |
endif |
611 |
|
612 |
! if (.not. haveRlist) then |
613 |
! write(default_error, *) 'rList has not been set in doForces!' |
614 |
! error = -1 |
615 |
! return |
616 |
! endif |
617 |
|
618 |
if (.not. haveNeighborList) then |
619 |
write(default_error, *) 'neighbor list has not been initialized in doForces!' |
620 |
error = -1 |
621 |
return |
622 |
end if |
623 |
|
624 |
if (.not. haveSaneForceField) then |
625 |
write(default_error, *) 'Force Field is not sane in doForces!' |
626 |
error = -1 |
627 |
return |
628 |
end if |
629 |
|
630 |
#ifdef IS_MPI |
631 |
if (.not. isMPISimSet()) then |
632 |
write(default_error,*) "ERROR: mpiSimulation has not been initialized!" |
633 |
error = -1 |
634 |
return |
635 |
endif |
636 |
#endif |
637 |
return |
638 |
end subroutine doReadyCheck |
639 |
|
640 |
|
641 |
subroutine init_FF(thisESM, thisStat) |
642 |
|
643 |
integer, intent(in) :: thisESM |
644 |
integer, intent(out) :: thisStat |
645 |
integer :: my_status, nMatches |
646 |
integer, pointer :: MatchList(:) => null() |
647 |
|
648 |
!! assume things are copacetic, unless they aren't |
649 |
thisStat = 0 |
650 |
|
651 |
electrostaticSummationMethod = thisESM |
652 |
|
653 |
!! init_FF is called *after* all of the atom types have been |
654 |
!! defined in atype_module using the new_atype subroutine. |
655 |
!! |
656 |
!! this will scan through the known atypes and figure out what |
657 |
!! interactions are used by the force field. |
658 |
|
659 |
FF_uses_DirectionalAtoms = .false. |
660 |
FF_uses_Dipoles = .false. |
661 |
FF_uses_GayBerne = .false. |
662 |
FF_uses_EAM = .false. |
663 |
|
664 |
call getMatchingElementList(atypes, "is_Directional", .true., & |
665 |
nMatches, MatchList) |
666 |
if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true. |
667 |
|
668 |
call getMatchingElementList(atypes, "is_Dipole", .true., & |
669 |
nMatches, MatchList) |
670 |
if (nMatches .gt. 0) FF_uses_Dipoles = .true. |
671 |
|
672 |
call getMatchingElementList(atypes, "is_GayBerne", .true., & |
673 |
nMatches, MatchList) |
674 |
if (nMatches .gt. 0) FF_uses_GayBerne = .true. |
675 |
|
676 |
call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList) |
677 |
if (nMatches .gt. 0) FF_uses_EAM = .true. |
678 |
|
679 |
|
680 |
haveSaneForceField = .true. |
681 |
|
682 |
if (FF_uses_EAM) then |
683 |
call init_EAM_FF(my_status) |
684 |
if (my_status /= 0) then |
685 |
write(default_error, *) "init_EAM_FF returned a bad status" |
686 |
thisStat = -1 |
687 |
haveSaneForceField = .false. |
688 |
return |
689 |
end if |
690 |
endif |
691 |
|
692 |
if (.not. haveNeighborList) then |
693 |
!! Create neighbor lists |
694 |
call expandNeighborList(nLocal, my_status) |
695 |
if (my_Status /= 0) then |
696 |
write(default_error,*) "SimSetup: ExpandNeighborList returned error." |
697 |
thisStat = -1 |
698 |
return |
699 |
endif |
700 |
haveNeighborList = .true. |
701 |
endif |
702 |
|
703 |
end subroutine init_FF |
704 |
|
705 |
|
706 |
!! Does force loop over i,j pairs. Calls do_pair to calculates forces. |
707 |
!-------------------------------------------------------------> |
708 |
subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, & |
709 |
do_pot_c, do_stress_c, error) |
710 |
!! Position array provided by C, dimensioned by getNlocal |
711 |
real ( kind = dp ), dimension(3, nLocal) :: q |
712 |
!! molecular center-of-mass position array |
713 |
real ( kind = dp ), dimension(3, nGroups) :: q_group |
714 |
!! Rotation Matrix for each long range particle in simulation. |
715 |
real( kind = dp), dimension(9, nLocal) :: A |
716 |
!! Unit vectors for dipoles (lab frame) |
717 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
718 |
!! Force array provided by C, dimensioned by getNlocal |
719 |
real ( kind = dp ), dimension(3,nLocal) :: f |
720 |
!! Torsion array provided by C, dimensioned by getNlocal |
721 |
real( kind = dp ), dimension(3,nLocal) :: t |
722 |
|
723 |
!! Stress Tensor |
724 |
real( kind = dp), dimension(9) :: tau |
725 |
real ( kind = dp ),dimension(LR_POT_TYPES) :: pot |
726 |
logical ( kind = 2) :: do_pot_c, do_stress_c |
727 |
logical :: do_pot |
728 |
logical :: do_stress |
729 |
logical :: in_switching_region |
730 |
#ifdef IS_MPI |
731 |
real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local |
732 |
integer :: nAtomsInRow |
733 |
integer :: nAtomsInCol |
734 |
integer :: nprocs |
735 |
integer :: nGroupsInRow |
736 |
integer :: nGroupsInCol |
737 |
#endif |
738 |
integer :: natoms |
739 |
logical :: update_nlist |
740 |
integer :: i, j, jstart, jend, jnab |
741 |
integer :: istart, iend |
742 |
integer :: ia, jb, atom1, atom2 |
743 |
integer :: nlist |
744 |
real( kind = DP ) :: ratmsq, rgrpsq, rgrp, vpair, vij |
745 |
real( kind = DP ) :: sw, dswdr, swderiv, mf |
746 |
real( kind = DP ) :: rVal |
747 |
real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij |
748 |
real(kind=dp), dimension(3) :: fstrs, f2strs |
749 |
real(kind=dp) :: rfpot, mu_i, virial |
750 |
integer :: me_i, me_j, n_in_i, n_in_j |
751 |
logical :: is_dp_i |
752 |
integer :: neighborListSize |
753 |
integer :: listerror, error |
754 |
integer :: localError |
755 |
integer :: propPack_i, propPack_j |
756 |
integer :: loopStart, loopEnd, loop |
757 |
integer :: iHash |
758 |
integer :: i1 |
759 |
|
760 |
|
761 |
!! initialize local variables |
762 |
|
763 |
#ifdef IS_MPI |
764 |
pot_local = 0.0_dp |
765 |
nAtomsInRow = getNatomsInRow(plan_atom_row) |
766 |
nAtomsInCol = getNatomsInCol(plan_atom_col) |
767 |
nGroupsInRow = getNgroupsInRow(plan_group_row) |
768 |
nGroupsInCol = getNgroupsInCol(plan_group_col) |
769 |
#else |
770 |
natoms = nlocal |
771 |
#endif |
772 |
|
773 |
call doReadyCheck(localError) |
774 |
if ( localError .ne. 0 ) then |
775 |
call handleError("do_force_loop", "Not Initialized") |
776 |
error = -1 |
777 |
return |
778 |
end if |
779 |
call zero_work_arrays() |
780 |
|
781 |
do_pot = do_pot_c |
782 |
do_stress = do_stress_c |
783 |
|
784 |
! Gather all information needed by all force loops: |
785 |
|
786 |
#ifdef IS_MPI |
787 |
|
788 |
call gather(q, q_Row, plan_atom_row_3d) |
789 |
call gather(q, q_Col, plan_atom_col_3d) |
790 |
|
791 |
call gather(q_group, q_group_Row, plan_group_row_3d) |
792 |
call gather(q_group, q_group_Col, plan_group_col_3d) |
793 |
|
794 |
if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then |
795 |
call gather(eFrame, eFrame_Row, plan_atom_row_rotation) |
796 |
call gather(eFrame, eFrame_Col, plan_atom_col_rotation) |
797 |
|
798 |
call gather(A, A_Row, plan_atom_row_rotation) |
799 |
call gather(A, A_Col, plan_atom_col_rotation) |
800 |
endif |
801 |
|
802 |
#endif |
803 |
|
804 |
!! Begin force loop timing: |
805 |
#ifdef PROFILE |
806 |
call cpu_time(forceTimeInitial) |
807 |
nloops = nloops + 1 |
808 |
#endif |
809 |
|
810 |
loopEnd = PAIR_LOOP |
811 |
if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then |
812 |
loopStart = PREPAIR_LOOP |
813 |
else |
814 |
loopStart = PAIR_LOOP |
815 |
endif |
816 |
|
817 |
do loop = loopStart, loopEnd |
818 |
|
819 |
! See if we need to update neighbor lists |
820 |
! (but only on the first time through): |
821 |
if (loop .eq. loopStart) then |
822 |
#ifdef IS_MPI |
823 |
call checkNeighborList(nGroupsInRow, q_group_row, listSkin, & |
824 |
update_nlist) |
825 |
#else |
826 |
call checkNeighborList(nGroups, q_group, listSkin, & |
827 |
update_nlist) |
828 |
#endif |
829 |
endif |
830 |
|
831 |
if (update_nlist) then |
832 |
!! save current configuration and construct neighbor list |
833 |
#ifdef IS_MPI |
834 |
call saveNeighborList(nGroupsInRow, q_group_row) |
835 |
#else |
836 |
call saveNeighborList(nGroups, q_group) |
837 |
#endif |
838 |
neighborListSize = size(list) |
839 |
nlist = 0 |
840 |
endif |
841 |
|
842 |
istart = 1 |
843 |
#ifdef IS_MPI |
844 |
iend = nGroupsInRow |
845 |
#else |
846 |
iend = nGroups - 1 |
847 |
#endif |
848 |
outer: do i = istart, iend |
849 |
|
850 |
if (update_nlist) point(i) = nlist + 1 |
851 |
|
852 |
n_in_i = groupStartRow(i+1) - groupStartRow(i) |
853 |
|
854 |
if (update_nlist) then |
855 |
#ifdef IS_MPI |
856 |
jstart = 1 |
857 |
jend = nGroupsInCol |
858 |
#else |
859 |
jstart = i+1 |
860 |
jend = nGroups |
861 |
#endif |
862 |
else |
863 |
jstart = point(i) |
864 |
jend = point(i+1) - 1 |
865 |
! make sure group i has neighbors |
866 |
if (jstart .gt. jend) cycle outer |
867 |
endif |
868 |
|
869 |
do jnab = jstart, jend |
870 |
if (update_nlist) then |
871 |
j = jnab |
872 |
else |
873 |
j = list(jnab) |
874 |
endif |
875 |
|
876 |
#ifdef IS_MPI |
877 |
me_j = atid_col(j) |
878 |
call get_interatomic_vector(q_group_Row(:,i), & |
879 |
q_group_Col(:,j), d_grp, rgrpsq) |
880 |
#else |
881 |
me_j = atid(j) |
882 |
call get_interatomic_vector(q_group(:,i), & |
883 |
q_group(:,j), d_grp, rgrpsq) |
884 |
#endif |
885 |
|
886 |
if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then |
887 |
if (update_nlist) then |
888 |
nlist = nlist + 1 |
889 |
|
890 |
if (nlist > neighborListSize) then |
891 |
#ifdef IS_MPI |
892 |
call expandNeighborList(nGroupsInRow, listerror) |
893 |
#else |
894 |
call expandNeighborList(nGroups, listerror) |
895 |
#endif |
896 |
if (listerror /= 0) then |
897 |
error = -1 |
898 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
899 |
return |
900 |
end if |
901 |
neighborListSize = size(list) |
902 |
endif |
903 |
|
904 |
list(nlist) = j |
905 |
endif |
906 |
|
907 |
if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then |
908 |
|
909 |
if (loop .eq. PAIR_LOOP) then |
910 |
vij = 0.0d0 |
911 |
fij(1:3) = 0.0d0 |
912 |
endif |
913 |
|
914 |
call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, & |
915 |
in_switching_region) |
916 |
|
917 |
n_in_j = groupStartCol(j+1) - groupStartCol(j) |
918 |
|
919 |
do ia = groupStartRow(i), groupStartRow(i+1)-1 |
920 |
|
921 |
atom1 = groupListRow(ia) |
922 |
|
923 |
inner: do jb = groupStartCol(j), groupStartCol(j+1)-1 |
924 |
|
925 |
atom2 = groupListCol(jb) |
926 |
|
927 |
if (skipThisPair(atom1, atom2)) cycle inner |
928 |
|
929 |
if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then |
930 |
d_atm(1:3) = d_grp(1:3) |
931 |
ratmsq = rgrpsq |
932 |
else |
933 |
#ifdef IS_MPI |
934 |
call get_interatomic_vector(q_Row(:,atom1), & |
935 |
q_Col(:,atom2), d_atm, ratmsq) |
936 |
#else |
937 |
call get_interatomic_vector(q(:,atom1), & |
938 |
q(:,atom2), d_atm, ratmsq) |
939 |
#endif |
940 |
endif |
941 |
|
942 |
if (loop .eq. PREPAIR_LOOP) then |
943 |
#ifdef IS_MPI |
944 |
call do_prepair(atom1, atom2, ratmsq, d_atm, sw, & |
945 |
rgrpsq, d_grp, do_pot, do_stress, & |
946 |
eFrame, A, f, t, pot_local) |
947 |
#else |
948 |
call do_prepair(atom1, atom2, ratmsq, d_atm, sw, & |
949 |
rgrpsq, d_grp, do_pot, do_stress, & |
950 |
eFrame, A, f, t, pot) |
951 |
#endif |
952 |
else |
953 |
#ifdef IS_MPI |
954 |
call do_pair(atom1, atom2, ratmsq, d_atm, sw, & |
955 |
do_pot, eFrame, A, f, t, pot_local, vpair, & |
956 |
fpair, d_grp, rgrp, fstrs) |
957 |
#else |
958 |
call do_pair(atom1, atom2, ratmsq, d_atm, sw, & |
959 |
do_pot, eFrame, A, f, t, pot, vpair, fpair, & |
960 |
d_grp, rgrp, fstrs) |
961 |
#endif |
962 |
f2strs(1:3) = f2strs(1:3) + fstrs(1:3) |
963 |
vij = vij + vpair |
964 |
fij(1:3) = fij(1:3) + fpair(1:3) |
965 |
endif |
966 |
enddo inner |
967 |
enddo |
968 |
|
969 |
if (loop .eq. PAIR_LOOP) then |
970 |
if (in_switching_region) then |
971 |
swderiv = vij*dswdr/rgrp |
972 |
fij(1) = fij(1) + swderiv*d_grp(1) |
973 |
fij(2) = fij(2) + swderiv*d_grp(2) |
974 |
fij(3) = fij(3) + swderiv*d_grp(3) |
975 |
|
976 |
do ia=groupStartRow(i), groupStartRow(i+1)-1 |
977 |
atom1=groupListRow(ia) |
978 |
mf = mfactRow(atom1) |
979 |
#ifdef IS_MPI |
980 |
f_Row(1,atom1) = f_Row(1,atom1) + swderiv*d_grp(1)*mf |
981 |
f_Row(2,atom1) = f_Row(2,atom1) + swderiv*d_grp(2)*mf |
982 |
f_Row(3,atom1) = f_Row(3,atom1) + swderiv*d_grp(3)*mf |
983 |
#else |
984 |
f(1,atom1) = f(1,atom1) + swderiv*d_grp(1)*mf |
985 |
f(2,atom1) = f(2,atom1) + swderiv*d_grp(2)*mf |
986 |
f(3,atom1) = f(3,atom1) + swderiv*d_grp(3)*mf |
987 |
#endif |
988 |
enddo |
989 |
|
990 |
do jb=groupStartCol(j), groupStartCol(j+1)-1 |
991 |
atom2=groupListCol(jb) |
992 |
mf = mfactCol(atom2) |
993 |
#ifdef IS_MPI |
994 |
f_Col(1,atom2) = f_Col(1,atom2) - swderiv*d_grp(1)*mf |
995 |
f_Col(2,atom2) = f_Col(2,atom2) - swderiv*d_grp(2)*mf |
996 |
f_Col(3,atom2) = f_Col(3,atom2) - swderiv*d_grp(3)*mf |
997 |
#else |
998 |
f(1,atom2) = f(1,atom2) - swderiv*d_grp(1)*mf |
999 |
f(2,atom2) = f(2,atom2) - swderiv*d_grp(2)*mf |
1000 |
f(3,atom2) = f(3,atom2) - swderiv*d_grp(3)*mf |
1001 |
#endif |
1002 |
enddo |
1003 |
endif |
1004 |
|
1005 |
if (do_stress) call add_stress_tensor(d_grp, fij) |
1006 |
endif |
1007 |
endif |
1008 |
endif |
1009 |
enddo |
1010 |
|
1011 |
enddo outer |
1012 |
|
1013 |
if (update_nlist) then |
1014 |
#ifdef IS_MPI |
1015 |
point(nGroupsInRow + 1) = nlist + 1 |
1016 |
#else |
1017 |
point(nGroups) = nlist + 1 |
1018 |
#endif |
1019 |
if (loop .eq. PREPAIR_LOOP) then |
1020 |
! we just did the neighbor list update on the first |
1021 |
! pass, so we don't need to do it |
1022 |
! again on the second pass |
1023 |
update_nlist = .false. |
1024 |
endif |
1025 |
endif |
1026 |
|
1027 |
if (loop .eq. PREPAIR_LOOP) then |
1028 |
call do_preforce(nlocal, pot) |
1029 |
endif |
1030 |
|
1031 |
enddo |
1032 |
|
1033 |
!! Do timing |
1034 |
#ifdef PROFILE |
1035 |
call cpu_time(forceTimeFinal) |
1036 |
forceTime = forceTime + forceTimeFinal - forceTimeInitial |
1037 |
#endif |
1038 |
|
1039 |
#ifdef IS_MPI |
1040 |
!!distribute forces |
1041 |
|
1042 |
f_temp = 0.0_dp |
1043 |
call scatter(f_Row,f_temp,plan_atom_row_3d) |
1044 |
do i = 1,nlocal |
1045 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
1046 |
end do |
1047 |
|
1048 |
f_temp = 0.0_dp |
1049 |
call scatter(f_Col,f_temp,plan_atom_col_3d) |
1050 |
do i = 1,nlocal |
1051 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
1052 |
end do |
1053 |
|
1054 |
if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then |
1055 |
t_temp = 0.0_dp |
1056 |
call scatter(t_Row,t_temp,plan_atom_row_3d) |
1057 |
do i = 1,nlocal |
1058 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
1059 |
end do |
1060 |
t_temp = 0.0_dp |
1061 |
call scatter(t_Col,t_temp,plan_atom_col_3d) |
1062 |
|
1063 |
do i = 1,nlocal |
1064 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
1065 |
end do |
1066 |
endif |
1067 |
|
1068 |
if (do_pot) then |
1069 |
! scatter/gather pot_row into the members of my column |
1070 |
do i = 1,LR_POT_TYPES |
1071 |
call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row) |
1072 |
end do |
1073 |
! scatter/gather pot_local into all other procs |
1074 |
! add resultant to get total pot |
1075 |
do i = 1, nlocal |
1076 |
pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) & |
1077 |
+ pot_Temp(1:LR_POT_TYPES,i) |
1078 |
enddo |
1079 |
|
1080 |
pot_Temp = 0.0_DP |
1081 |
do i = 1,LR_POT_TYPES |
1082 |
call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col) |
1083 |
end do |
1084 |
do i = 1, nlocal |
1085 |
pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)& |
1086 |
+ pot_Temp(1:LR_POT_TYPES,i) |
1087 |
enddo |
1088 |
|
1089 |
endif |
1090 |
#endif |
1091 |
|
1092 |
if (SIM_requires_postpair_calc) then |
1093 |
do i = 1, nlocal |
1094 |
|
1095 |
! we loop only over the local atoms, so we don't need row and column |
1096 |
! lookups for the types |
1097 |
|
1098 |
me_i = atid(i) |
1099 |
|
1100 |
! is the atom electrostatic? See if it would have an |
1101 |
! electrostatic interaction with itself |
1102 |
iHash = InteractionHash(me_i,me_i) |
1103 |
|
1104 |
if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then |
1105 |
#ifdef IS_MPI |
1106 |
call self_self(i, eFrame, pot_local(ELECTROSTATIC_POT), & |
1107 |
t, do_pot) |
1108 |
#else |
1109 |
call self_self(i, eFrame, pot(ELECTROSTATIC_POT), & |
1110 |
t, do_pot) |
1111 |
#endif |
1112 |
endif |
1113 |
|
1114 |
|
1115 |
if (electrostaticSummationMethod.eq.REACTION_FIELD) then |
1116 |
|
1117 |
! loop over the excludes to accumulate RF stuff we've |
1118 |
! left out of the normal pair loop |
1119 |
|
1120 |
do i1 = 1, nSkipsForAtom(i) |
1121 |
j = skipsForAtom(i, i1) |
1122 |
|
1123 |
! prevent overcounting of the skips |
1124 |
if (i.lt.j) then |
1125 |
call get_interatomic_vector(q(:,i), & |
1126 |
q(:,j), d_atm, ratmsq) |
1127 |
rVal = dsqrt(ratmsq) |
1128 |
call get_switch(ratmsq, sw, dswdr, rVal, group_switch, & |
1129 |
in_switching_region) |
1130 |
#ifdef IS_MPI |
1131 |
call rf_self_excludes(i, j, sw, eFrame, d_atm, rVal, & |
1132 |
vpair, pot_local(ELECTROSTATIC_POT), f, t, do_pot) |
1133 |
#else |
1134 |
call rf_self_excludes(i, j, sw, eFrame, d_atm, rVal, & |
1135 |
vpair, pot(ELECTROSTATIC_POT), f, t, do_pot) |
1136 |
#endif |
1137 |
endif |
1138 |
enddo |
1139 |
endif |
1140 |
enddo |
1141 |
endif |
1142 |
|
1143 |
#ifdef IS_MPI |
1144 |
|
1145 |
if (do_pot) then |
1146 |
call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision,mpi_sum, & |
1147 |
mpi_comm_world,mpi_err) |
1148 |
endif |
1149 |
|
1150 |
if (do_stress) then |
1151 |
call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, & |
1152 |
mpi_comm_world,mpi_err) |
1153 |
call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, & |
1154 |
mpi_comm_world,mpi_err) |
1155 |
endif |
1156 |
|
1157 |
#else |
1158 |
|
1159 |
if (do_stress) then |
1160 |
tau = tau_Temp |
1161 |
virial = virial_Temp |
1162 |
endif |
1163 |
|
1164 |
#endif |
1165 |
|
1166 |
end subroutine do_force_loop |
1167 |
|
1168 |
!!$ subroutine do_pair(i, j, rijsq, d, sw, do_pot, & |
1169 |
!!$ eFrame, A, f, t, pot, vpair, fpair, d_grp, r_grp) |
1170 |
subroutine do_pair(i, j, rijsq, d, sw, do_pot, & |
1171 |
eFrame, A, f, t, pot, vpair, fpair, d_grp, r_grp, fstrs) |
1172 |
|
1173 |
real( kind = dp ) :: vpair, sw |
1174 |
real( kind = dp ), dimension(LR_POT_TYPES) :: pot |
1175 |
real( kind = dp ), dimension(3) :: fpair |
1176 |
real( kind = dp ), dimension(3) :: fstrs |
1177 |
real( kind = dp ), dimension(nLocal) :: mfact |
1178 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
1179 |
real( kind = dp ), dimension(9,nLocal) :: A |
1180 |
real( kind = dp ), dimension(3,nLocal) :: f |
1181 |
real( kind = dp ), dimension(3,nLocal) :: t |
1182 |
|
1183 |
logical, intent(inout) :: do_pot |
1184 |
integer, intent(in) :: i, j |
1185 |
real ( kind = dp ), intent(inout) :: rijsq |
1186 |
real ( kind = dp ), intent(inout) :: r_grp |
1187 |
real ( kind = dp ), intent(inout) :: d(3) |
1188 |
real ( kind = dp ), intent(inout) :: d_grp(3) |
1189 |
real ( kind = dp ) :: r |
1190 |
integer :: me_i, me_j |
1191 |
|
1192 |
integer :: iHash |
1193 |
|
1194 |
r = sqrt(rijsq) |
1195 |
vpair = 0.0d0 |
1196 |
fpair(1:3) = 0.0d0 |
1197 |
|
1198 |
#ifdef IS_MPI |
1199 |
me_i = atid_row(i) |
1200 |
me_j = atid_col(j) |
1201 |
#else |
1202 |
me_i = atid(i) |
1203 |
me_j = atid(j) |
1204 |
#endif |
1205 |
|
1206 |
iHash = InteractionHash(me_i, me_j) |
1207 |
|
1208 |
if ( iand(iHash, LJ_PAIR).ne.0 ) then |
1209 |
call do_lj_pair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1210 |
pot(VDW_POT), f, do_pot) |
1211 |
endif |
1212 |
|
1213 |
if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then |
1214 |
!!$ call doElectrostaticPair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1215 |
!!$ pot(ELECTROSTATIC_POT), eFrame, f, t, do_pot) |
1216 |
call doElectrostaticPair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1217 |
pot(ELECTROSTATIC_POT), eFrame, f, t, do_pot, fstrs) |
1218 |
endif |
1219 |
|
1220 |
if ( iand(iHash, STICKY_PAIR).ne.0 ) then |
1221 |
call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1222 |
pot(HB_POT), A, f, t, do_pot) |
1223 |
endif |
1224 |
|
1225 |
if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then |
1226 |
call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1227 |
pot(HB_POT), A, f, t, do_pot) |
1228 |
endif |
1229 |
|
1230 |
if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then |
1231 |
call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1232 |
pot(VDW_POT), A, f, t, do_pot) |
1233 |
endif |
1234 |
|
1235 |
if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then |
1236 |
call do_gb_lj_pair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1237 |
pot(VDW_POT), A, f, t, do_pot) |
1238 |
endif |
1239 |
|
1240 |
if ( iand(iHash, EAM_PAIR).ne.0 ) then |
1241 |
call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1242 |
pot(METALLIC_POT), f, do_pot) |
1243 |
endif |
1244 |
|
1245 |
if ( iand(iHash, SHAPE_PAIR).ne.0 ) then |
1246 |
call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1247 |
pot(VDW_POT), A, f, t, do_pot) |
1248 |
endif |
1249 |
|
1250 |
if ( iand(iHash, SHAPE_LJ).ne.0 ) then |
1251 |
call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, & |
1252 |
pot(VDW_POT), A, f, t, do_pot) |
1253 |
endif |
1254 |
|
1255 |
end subroutine do_pair |
1256 |
|
1257 |
subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, & |
1258 |
do_pot, do_stress, eFrame, A, f, t, pot) |
1259 |
|
1260 |
real( kind = dp ) :: sw |
1261 |
real( kind = dp ), dimension(LR_POT_TYPES) :: pot |
1262 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
1263 |
real (kind=dp), dimension(9,nLocal) :: A |
1264 |
real (kind=dp), dimension(3,nLocal) :: f |
1265 |
real (kind=dp), dimension(3,nLocal) :: t |
1266 |
|
1267 |
logical, intent(inout) :: do_pot, do_stress |
1268 |
integer, intent(in) :: i, j |
1269 |
real ( kind = dp ), intent(inout) :: rijsq, rcijsq |
1270 |
real ( kind = dp ) :: r, rc |
1271 |
real ( kind = dp ), intent(inout) :: d(3), dc(3) |
1272 |
|
1273 |
integer :: me_i, me_j, iHash |
1274 |
|
1275 |
r = sqrt(rijsq) |
1276 |
|
1277 |
#ifdef IS_MPI |
1278 |
me_i = atid_row(i) |
1279 |
me_j = atid_col(j) |
1280 |
#else |
1281 |
me_i = atid(i) |
1282 |
me_j = atid(j) |
1283 |
#endif |
1284 |
|
1285 |
iHash = InteractionHash(me_i, me_j) |
1286 |
|
1287 |
if ( iand(iHash, EAM_PAIR).ne.0 ) then |
1288 |
call calc_EAM_prepair_rho(i, j, d, r, rijsq ) |
1289 |
endif |
1290 |
|
1291 |
end subroutine do_prepair |
1292 |
|
1293 |
|
1294 |
subroutine do_preforce(nlocal,pot) |
1295 |
integer :: nlocal |
1296 |
real( kind = dp ),dimension(LR_POT_TYPES) :: pot |
1297 |
|
1298 |
if (FF_uses_EAM .and. SIM_uses_EAM) then |
1299 |
call calc_EAM_preforce_Frho(nlocal,pot(METALLIC_POT)) |
1300 |
endif |
1301 |
|
1302 |
|
1303 |
end subroutine do_preforce |
1304 |
|
1305 |
|
1306 |
subroutine get_interatomic_vector(q_i, q_j, d, r_sq) |
1307 |
|
1308 |
real (kind = dp), dimension(3) :: q_i |
1309 |
real (kind = dp), dimension(3) :: q_j |
1310 |
real ( kind = dp ), intent(out) :: r_sq |
1311 |
real( kind = dp ) :: d(3), scaled(3) |
1312 |
integer i |
1313 |
|
1314 |
d(1:3) = q_j(1:3) - q_i(1:3) |
1315 |
|
1316 |
! Wrap back into periodic box if necessary |
1317 |
if ( SIM_uses_PBC ) then |
1318 |
|
1319 |
if( .not.boxIsOrthorhombic ) then |
1320 |
! calc the scaled coordinates. |
1321 |
|
1322 |
scaled = matmul(HmatInv, d) |
1323 |
|
1324 |
! wrap the scaled coordinates |
1325 |
|
1326 |
scaled = scaled - anint(scaled) |
1327 |
|
1328 |
|
1329 |
! calc the wrapped real coordinates from the wrapped scaled |
1330 |
! coordinates |
1331 |
|
1332 |
d = matmul(Hmat,scaled) |
1333 |
|
1334 |
else |
1335 |
! calc the scaled coordinates. |
1336 |
|
1337 |
do i = 1, 3 |
1338 |
scaled(i) = d(i) * HmatInv(i,i) |
1339 |
|
1340 |
! wrap the scaled coordinates |
1341 |
|
1342 |
scaled(i) = scaled(i) - anint(scaled(i)) |
1343 |
|
1344 |
! calc the wrapped real coordinates from the wrapped scaled |
1345 |
! coordinates |
1346 |
|
1347 |
d(i) = scaled(i)*Hmat(i,i) |
1348 |
enddo |
1349 |
endif |
1350 |
|
1351 |
endif |
1352 |
|
1353 |
r_sq = dot_product(d,d) |
1354 |
|
1355 |
end subroutine get_interatomic_vector |
1356 |
|
1357 |
subroutine zero_work_arrays() |
1358 |
|
1359 |
#ifdef IS_MPI |
1360 |
|
1361 |
q_Row = 0.0_dp |
1362 |
q_Col = 0.0_dp |
1363 |
|
1364 |
q_group_Row = 0.0_dp |
1365 |
q_group_Col = 0.0_dp |
1366 |
|
1367 |
eFrame_Row = 0.0_dp |
1368 |
eFrame_Col = 0.0_dp |
1369 |
|
1370 |
A_Row = 0.0_dp |
1371 |
A_Col = 0.0_dp |
1372 |
|
1373 |
f_Row = 0.0_dp |
1374 |
f_Col = 0.0_dp |
1375 |
f_Temp = 0.0_dp |
1376 |
|
1377 |
t_Row = 0.0_dp |
1378 |
t_Col = 0.0_dp |
1379 |
t_Temp = 0.0_dp |
1380 |
|
1381 |
pot_Row = 0.0_dp |
1382 |
pot_Col = 0.0_dp |
1383 |
pot_Temp = 0.0_dp |
1384 |
|
1385 |
#endif |
1386 |
|
1387 |
if (FF_uses_EAM .and. SIM_uses_EAM) then |
1388 |
call clean_EAM() |
1389 |
endif |
1390 |
|
1391 |
tau_Temp = 0.0_dp |
1392 |
virial_Temp = 0.0_dp |
1393 |
end subroutine zero_work_arrays |
1394 |
|
1395 |
function skipThisPair(atom1, atom2) result(skip_it) |
1396 |
integer, intent(in) :: atom1 |
1397 |
integer, intent(in), optional :: atom2 |
1398 |
logical :: skip_it |
1399 |
integer :: unique_id_1, unique_id_2 |
1400 |
integer :: me_i,me_j |
1401 |
integer :: i |
1402 |
|
1403 |
skip_it = .false. |
1404 |
|
1405 |
!! there are a number of reasons to skip a pair or a particle |
1406 |
!! mostly we do this to exclude atoms who are involved in short |
1407 |
!! range interactions (bonds, bends, torsions), but we also need |
1408 |
!! to exclude some overcounted interactions that result from |
1409 |
!! the parallel decomposition |
1410 |
|
1411 |
#ifdef IS_MPI |
1412 |
!! in MPI, we have to look up the unique IDs for each atom |
1413 |
unique_id_1 = AtomRowToGlobal(atom1) |
1414 |
#else |
1415 |
!! in the normal loop, the atom numbers are unique |
1416 |
unique_id_1 = atom1 |
1417 |
#endif |
1418 |
|
1419 |
!! We were called with only one atom, so just check the global exclude |
1420 |
!! list for this atom |
1421 |
if (.not. present(atom2)) then |
1422 |
do i = 1, nExcludes_global |
1423 |
if (excludesGlobal(i) == unique_id_1) then |
1424 |
skip_it = .true. |
1425 |
return |
1426 |
end if |
1427 |
end do |
1428 |
return |
1429 |
end if |
1430 |
|
1431 |
#ifdef IS_MPI |
1432 |
unique_id_2 = AtomColToGlobal(atom2) |
1433 |
#else |
1434 |
unique_id_2 = atom2 |
1435 |
#endif |
1436 |
|
1437 |
#ifdef IS_MPI |
1438 |
!! this situation should only arise in MPI simulations |
1439 |
if (unique_id_1 == unique_id_2) then |
1440 |
skip_it = .true. |
1441 |
return |
1442 |
end if |
1443 |
|
1444 |
!! this prevents us from doing the pair on multiple processors |
1445 |
if (unique_id_1 < unique_id_2) then |
1446 |
if (mod(unique_id_1 + unique_id_2,2) == 0) then |
1447 |
skip_it = .true. |
1448 |
return |
1449 |
endif |
1450 |
else |
1451 |
if (mod(unique_id_1 + unique_id_2,2) == 1) then |
1452 |
skip_it = .true. |
1453 |
return |
1454 |
endif |
1455 |
endif |
1456 |
#endif |
1457 |
|
1458 |
!! the rest of these situations can happen in all simulations: |
1459 |
do i = 1, nExcludes_global |
1460 |
if ((excludesGlobal(i) == unique_id_1) .or. & |
1461 |
(excludesGlobal(i) == unique_id_2)) then |
1462 |
skip_it = .true. |
1463 |
return |
1464 |
endif |
1465 |
enddo |
1466 |
|
1467 |
do i = 1, nSkipsForAtom(atom1) |
1468 |
if (skipsForAtom(atom1, i) .eq. unique_id_2) then |
1469 |
skip_it = .true. |
1470 |
return |
1471 |
endif |
1472 |
end do |
1473 |
|
1474 |
return |
1475 |
end function skipThisPair |
1476 |
|
1477 |
function FF_UsesDirectionalAtoms() result(doesit) |
1478 |
logical :: doesit |
1479 |
doesit = FF_uses_DirectionalAtoms |
1480 |
end function FF_UsesDirectionalAtoms |
1481 |
|
1482 |
function FF_RequiresPrepairCalc() result(doesit) |
1483 |
logical :: doesit |
1484 |
doesit = FF_uses_EAM |
1485 |
end function FF_RequiresPrepairCalc |
1486 |
|
1487 |
#ifdef PROFILE |
1488 |
function getforcetime() result(totalforcetime) |
1489 |
real(kind=dp) :: totalforcetime |
1490 |
totalforcetime = forcetime |
1491 |
end function getforcetime |
1492 |
#endif |
1493 |
|
1494 |
!! This cleans componets of force arrays belonging only to fortran |
1495 |
|
1496 |
subroutine add_stress_tensor(dpair, fpair) |
1497 |
|
1498 |
real( kind = dp ), dimension(3), intent(in) :: dpair, fpair |
1499 |
|
1500 |
! because the d vector is the rj - ri vector, and |
1501 |
! because fx, fy, fz are the force on atom i, we need a |
1502 |
! negative sign here: |
1503 |
|
1504 |
tau_Temp(1) = tau_Temp(1) - dpair(1) * fpair(1) |
1505 |
tau_Temp(2) = tau_Temp(2) - dpair(1) * fpair(2) |
1506 |
tau_Temp(3) = tau_Temp(3) - dpair(1) * fpair(3) |
1507 |
tau_Temp(4) = tau_Temp(4) - dpair(2) * fpair(1) |
1508 |
tau_Temp(5) = tau_Temp(5) - dpair(2) * fpair(2) |
1509 |
tau_Temp(6) = tau_Temp(6) - dpair(2) * fpair(3) |
1510 |
tau_Temp(7) = tau_Temp(7) - dpair(3) * fpair(1) |
1511 |
tau_Temp(8) = tau_Temp(8) - dpair(3) * fpair(2) |
1512 |
tau_Temp(9) = tau_Temp(9) - dpair(3) * fpair(3) |
1513 |
|
1514 |
virial_Temp = virial_Temp + & |
1515 |
(tau_Temp(1) + tau_Temp(5) + tau_Temp(9)) |
1516 |
|
1517 |
end subroutine add_stress_tensor |
1518 |
|
1519 |
end module doForces |