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