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