| 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 status |
| 61 |
use ISO_C_BINDING |
| 62 |
|
| 63 |
#ifdef IS_MPI |
| 64 |
use mpiSimulation |
| 65 |
#endif |
| 66 |
|
| 67 |
implicit none |
| 68 |
PRIVATE |
| 69 |
|
| 70 |
real(kind=dp), external :: get_cutoff |
| 71 |
|
| 72 |
#define __FORTRAN90 |
| 73 |
#include "UseTheForce/DarkSide/fInteractionMap.h" |
| 74 |
#include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h" |
| 75 |
|
| 76 |
INTEGER, PARAMETER:: PREPAIR_LOOP = 1 |
| 77 |
INTEGER, PARAMETER:: PAIR_LOOP = 2 |
| 78 |
|
| 79 |
logical, save :: haveNeighborList = .false. |
| 80 |
logical, save :: haveSIMvariables = .false. |
| 81 |
logical, save :: haveSaneForceField = .false. |
| 82 |
logical, save :: haveDefaultCutoffs = .false. |
| 83 |
logical, save :: haveSkinThickness = .false. |
| 84 |
logical, save :: haveElectrostaticSummationMethod = .false. |
| 85 |
|
| 86 |
logical, save :: FF_uses_DirectionalAtoms |
| 87 |
logical, save :: FF_uses_Dipoles |
| 88 |
logical, save :: FF_uses_GayBerne |
| 89 |
logical, save :: FF_uses_EAM |
| 90 |
logical, save :: FF_uses_SC |
| 91 |
logical, save :: FF_uses_MNM |
| 92 |
|
| 93 |
logical, save :: SIM_uses_DirectionalAtoms |
| 94 |
logical, save :: SIM_uses_EAM |
| 95 |
logical, save :: SIM_uses_SC |
| 96 |
logical, save :: SIM_uses_MNM |
| 97 |
logical, save :: SIM_requires_postpair_calc |
| 98 |
logical, save :: SIM_requires_prepair_calc |
| 99 |
logical, save :: SIM_uses_PBC |
| 100 |
logical, save :: SIM_uses_AtomicVirial |
| 101 |
|
| 102 |
integer, save :: electrostaticSummationMethod |
| 103 |
|
| 104 |
real(kind=dp), save :: defaultRcut, defaultRsw, largestRcut |
| 105 |
real(kind=dp), save :: skinThickness |
| 106 |
logical, save :: defaultDoShiftPot |
| 107 |
logical, save :: defaultDoShiftFrc |
| 108 |
|
| 109 |
public :: init_FF |
| 110 |
public :: setCutoffs |
| 111 |
public :: setElectrostaticMethod |
| 112 |
public :: setSkinThickness |
| 113 |
public :: do_force_loop |
| 114 |
|
| 115 |
#ifdef PROFILE |
| 116 |
public :: getforcetime |
| 117 |
real, save :: forceTime = 0 |
| 118 |
real :: forceTimeInitial, forceTimeFinal |
| 119 |
integer :: nLoops |
| 120 |
#endif |
| 121 |
|
| 122 |
contains |
| 123 |
|
| 124 |
subroutine setCutoffs(defRcut, defRsw, defSP, defSF) |
| 125 |
|
| 126 |
real(kind=dp),intent(in) :: defRcut, defRsw |
| 127 |
integer, intent(in) :: defSP, defSF |
| 128 |
character(len = statusMsgSize) :: errMsg |
| 129 |
integer :: localError |
| 130 |
|
| 131 |
defaultRcut = defRcut |
| 132 |
defaultRsw = defRsw |
| 133 |
|
| 134 |
if (defSP .ne. 0) then |
| 135 |
defaultDoShiftPot = .true. |
| 136 |
else |
| 137 |
defaultDoShiftPot = .false. |
| 138 |
endif |
| 139 |
if (defSF .ne. 0) then |
| 140 |
defaultDoShiftFrc = .true. |
| 141 |
else |
| 142 |
defaultDoShiftFrc = .false. |
| 143 |
endif |
| 144 |
|
| 145 |
if (abs(defaultRcut-defaultRsw) .lt. 0.0001) then |
| 146 |
if (defaultDoShiftFrc) then |
| 147 |
write(errMsg, *) & |
| 148 |
'cutoffRadius and switchingRadius are set to the', newline & |
| 149 |
// tab, 'same value. OpenMD will use shifted force', newline & |
| 150 |
// tab, 'potentials instead of switching functions.' |
| 151 |
|
| 152 |
call handleInfo("setCutoffs", errMsg) |
| 153 |
else |
| 154 |
write(errMsg, *) & |
| 155 |
'cutoffRadius and switchingRadius are set to the', newline & |
| 156 |
// tab, 'same value. OpenMD will use shifted', newline & |
| 157 |
// tab, 'potentials instead of switching functions.' |
| 158 |
|
| 159 |
call handleInfo("setCutoffs", errMsg) |
| 160 |
|
| 161 |
defaultDoShiftPot = .true. |
| 162 |
endif |
| 163 |
|
| 164 |
endif |
| 165 |
|
| 166 |
localError = 0 |
| 167 |
call set_switch(defaultRsw, defaultRcut) |
| 168 |
call setHmatDangerousRcutValue(defaultRcut) |
| 169 |
|
| 170 |
haveDefaultCutoffs = .true. |
| 171 |
end subroutine setCutoffs |
| 172 |
|
| 173 |
subroutine setElectrostaticMethod( thisESM ) |
| 174 |
|
| 175 |
integer, intent(in) :: thisESM |
| 176 |
|
| 177 |
electrostaticSummationMethod = thisESM |
| 178 |
haveElectrostaticSummationMethod = .true. |
| 179 |
|
| 180 |
end subroutine setElectrostaticMethod |
| 181 |
|
| 182 |
subroutine setSkinThickness( thisSkin ) |
| 183 |
|
| 184 |
real(kind=dp), intent(in) :: thisSkin |
| 185 |
|
| 186 |
skinThickness = thisSkin |
| 187 |
haveSkinThickness = .true. |
| 188 |
haveGtypeCutoffMap = .false. |
| 189 |
|
| 190 |
end subroutine setSkinThickness |
| 191 |
|
| 192 |
subroutine setSimVariables() |
| 193 |
SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms() |
| 194 |
SIM_uses_EAM = SimUsesEAM() |
| 195 |
SIM_requires_postpair_calc = SimRequiresPostpairCalc() |
| 196 |
SIM_requires_prepair_calc = SimRequiresPrepairCalc() |
| 197 |
SIM_uses_PBC = SimUsesPBC() |
| 198 |
SIM_uses_SC = SimUsesSC() |
| 199 |
SIM_uses_AtomicVirial = SimUsesAtomicVirial() |
| 200 |
|
| 201 |
haveSIMvariables = .true. |
| 202 |
|
| 203 |
return |
| 204 |
end subroutine setSimVariables |
| 205 |
|
| 206 |
subroutine doReadyCheck(error) |
| 207 |
integer, intent(out) :: error |
| 208 |
integer :: myStatus |
| 209 |
|
| 210 |
error = 0 |
| 211 |
|
| 212 |
if (.not. haveSIMvariables) then |
| 213 |
call setSimVariables() |
| 214 |
endif |
| 215 |
|
| 216 |
if (.not. haveNeighborList) then |
| 217 |
write(default_error, *) 'neighbor list has not been initialized in doForces!' |
| 218 |
error = -1 |
| 219 |
return |
| 220 |
end if |
| 221 |
|
| 222 |
if (.not. haveSaneForceField) then |
| 223 |
write(default_error, *) 'Force Field is not sane in doForces!' |
| 224 |
error = -1 |
| 225 |
return |
| 226 |
end if |
| 227 |
|
| 228 |
#ifdef IS_MPI |
| 229 |
if (.not. isMPISimSet()) then |
| 230 |
write(default_error,*) "ERROR: mpiSimulation has not been initialized!" |
| 231 |
error = -1 |
| 232 |
return |
| 233 |
endif |
| 234 |
#endif |
| 235 |
return |
| 236 |
end subroutine doReadyCheck |
| 237 |
|
| 238 |
subroutine init_FF(thisStat) |
| 239 |
|
| 240 |
integer, intent(out) :: thisStat |
| 241 |
integer :: my_status, nMatches |
| 242 |
integer, pointer :: MatchList(:) => null() |
| 243 |
|
| 244 |
!! assume things are copacetic, unless they aren't |
| 245 |
thisStat = 0 |
| 246 |
|
| 247 |
!! init_FF is called *after* all of the atom types have been |
| 248 |
!! defined in atype_module using the new_atype subroutine. |
| 249 |
!! |
| 250 |
!! this will scan through the known atypes and figure out what |
| 251 |
!! interactions are used by the force field. |
| 252 |
|
| 253 |
FF_uses_DirectionalAtoms = .false. |
| 254 |
FF_uses_Dipoles = .false. |
| 255 |
FF_uses_GayBerne = .false. |
| 256 |
FF_uses_EAM = .false. |
| 257 |
FF_uses_SC = .false. |
| 258 |
|
| 259 |
call getMatchingElementList(atypes, "is_Directional", .true., & |
| 260 |
nMatches, MatchList) |
| 261 |
if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true. |
| 262 |
|
| 263 |
call getMatchingElementList(atypes, "is_Dipole", .true., & |
| 264 |
nMatches, MatchList) |
| 265 |
if (nMatches .gt. 0) FF_uses_Dipoles = .true. |
| 266 |
|
| 267 |
call getMatchingElementList(atypes, "is_GayBerne", .true., & |
| 268 |
nMatches, MatchList) |
| 269 |
if (nMatches .gt. 0) FF_uses_GayBerne = .true. |
| 270 |
|
| 271 |
call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList) |
| 272 |
if (nMatches .gt. 0) FF_uses_EAM = .true. |
| 273 |
|
| 274 |
call getMatchingElementList(atypes, "is_SC", .true., nMatches, MatchList) |
| 275 |
if (nMatches .gt. 0) FF_uses_SC = .true. |
| 276 |
|
| 277 |
|
| 278 |
haveSaneForceField = .true. |
| 279 |
|
| 280 |
|
| 281 |
if (.not. haveNeighborList) then |
| 282 |
!! Create neighbor lists |
| 283 |
call expandNeighborList(nLocal, my_status) |
| 284 |
if (my_Status /= 0) then |
| 285 |
write(default_error,*) "SimSetup: ExpandNeighborList returned error." |
| 286 |
thisStat = -1 |
| 287 |
return |
| 288 |
endif |
| 289 |
haveNeighborList = .true. |
| 290 |
endif |
| 291 |
|
| 292 |
end subroutine init_FF |
| 293 |
|
| 294 |
|
| 295 |
!! Does force loop over i,j pairs. Calls do_pair to calculates forces. |
| 296 |
!-------------------------------------------------------------> |
| 297 |
subroutine do_force_loop(q, q_group, A, eFrame, f, t, tau, pot, particle_pot, & |
| 298 |
error) |
| 299 |
!! Position array provided by C, dimensioned by getNlocal |
| 300 |
real ( kind = dp ), dimension(3, nLocal) :: q |
| 301 |
!! molecular center-of-mass position array |
| 302 |
real ( kind = dp ), dimension(3, nGroups) :: q_group |
| 303 |
!! Rotation Matrix for each long range particle in simulation. |
| 304 |
real( kind = dp), dimension(9, nLocal) :: A |
| 305 |
!! Unit vectors for dipoles (lab frame) |
| 306 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
| 307 |
!! Force array provided by C, dimensioned by getNlocal |
| 308 |
real ( kind = dp ), dimension(3,nLocal) :: f |
| 309 |
!! Torsion array provided by C, dimensioned by getNlocal |
| 310 |
real( kind = dp ), dimension(3,nLocal) :: t |
| 311 |
|
| 312 |
!! Stress Tensor |
| 313 |
real( kind = dp), dimension(9) :: tau |
| 314 |
real ( kind = dp ),dimension(LR_POT_TYPES) :: pot |
| 315 |
real( kind = dp ), dimension(nLocal) :: particle_pot |
| 316 |
real( kind = dp ), dimension(nLocal) :: skipped_charge |
| 317 |
|
| 318 |
logical :: in_switching_region |
| 319 |
#ifdef IS_MPI |
| 320 |
real( kind = DP ), dimension(LR_POT_TYPES) :: pot_local |
| 321 |
integer :: nAtomsInRow |
| 322 |
integer :: nAtomsInCol |
| 323 |
integer :: nprocs |
| 324 |
integer :: nGroupsInRow |
| 325 |
integer :: nGroupsInCol |
| 326 |
#endif |
| 327 |
integer :: natoms |
| 328 |
logical :: update_nlist |
| 329 |
integer :: i, j, jstart, jend, jnab |
| 330 |
integer :: istart, iend |
| 331 |
integer :: ia, jb, atom1, atom2 |
| 332 |
integer :: nlist |
| 333 |
real( kind = DP ) :: ratmsq, rgrpsq, rgrp, rag, vpair, vij |
| 334 |
real( kind = DP ) :: sw, dswdr, swderiv, mf |
| 335 |
real( kind = DP ) :: rVal |
| 336 |
real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij, fg, dag |
| 337 |
real(kind=dp) :: rfpot, mu_i |
| 338 |
real(kind=dp):: rCut |
| 339 |
integer :: me_i, me_j, n_in_i, n_in_j, iG, j1 |
| 340 |
logical :: is_dp_i |
| 341 |
integer :: neighborListSize |
| 342 |
integer :: listerror, error |
| 343 |
integer :: localError |
| 344 |
integer :: propPack_i, propPack_j |
| 345 |
integer :: loopStart, loopEnd, loop |
| 346 |
integer :: i1, topoDist |
| 347 |
|
| 348 |
real(kind=dp) :: skch |
| 349 |
|
| 350 |
!! initialize local variables |
| 351 |
|
| 352 |
#ifdef IS_MPI |
| 353 |
pot_local = 0.0_dp |
| 354 |
nAtomsInRow = getNatomsInRow(plan_atom_row) |
| 355 |
nAtomsInCol = getNatomsInCol(plan_atom_col) |
| 356 |
nGroupsInRow = getNgroupsInRow(plan_group_row) |
| 357 |
nGroupsInCol = getNgroupsInCol(plan_group_col) |
| 358 |
#else |
| 359 |
natoms = nlocal |
| 360 |
#endif |
| 361 |
|
| 362 |
call doReadyCheck(localError) |
| 363 |
if ( localError .ne. 0 ) then |
| 364 |
call handleError("do_force_loop", "Not Initialized") |
| 365 |
error = -1 |
| 366 |
return |
| 367 |
end if |
| 368 |
call zero_work_arrays() |
| 369 |
|
| 370 |
! Gather all information needed by all force loops: |
| 371 |
|
| 372 |
#ifdef IS_MPI |
| 373 |
|
| 374 |
call gather(q, q_Row, plan_atom_row_3d) |
| 375 |
call gather(q, q_Col, plan_atom_col_3d) |
| 376 |
|
| 377 |
call gather(q_group, q_group_Row, plan_group_row_3d) |
| 378 |
call gather(q_group, q_group_Col, plan_group_col_3d) |
| 379 |
|
| 380 |
if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then |
| 381 |
call gather(eFrame, eFrame_Row, plan_atom_row_rotation) |
| 382 |
call gather(eFrame, eFrame_Col, plan_atom_col_rotation) |
| 383 |
|
| 384 |
call gather(A, A_Row, plan_atom_row_rotation) |
| 385 |
call gather(A, A_Col, plan_atom_col_rotation) |
| 386 |
endif |
| 387 |
|
| 388 |
#endif |
| 389 |
|
| 390 |
!! Begin force loop timing: |
| 391 |
#ifdef PROFILE |
| 392 |
call cpu_time(forceTimeInitial) |
| 393 |
nloops = nloops + 1 |
| 394 |
#endif |
| 395 |
|
| 396 |
loopEnd = PAIR_LOOP |
| 397 |
if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then |
| 398 |
loopStart = PREPAIR_LOOP |
| 399 |
else |
| 400 |
loopStart = PAIR_LOOP |
| 401 |
endif |
| 402 |
|
| 403 |
do loop = loopStart, loopEnd |
| 404 |
|
| 405 |
! See if we need to update neighbor lists |
| 406 |
! (but only on the first time through): |
| 407 |
if (loop .eq. loopStart) then |
| 408 |
#ifdef IS_MPI |
| 409 |
call checkNeighborList(nGroupsInRow, q_group_row, skinThickness, & |
| 410 |
update_nlist) |
| 411 |
#else |
| 412 |
call checkNeighborList(nGroups, q_group, skinThickness, & |
| 413 |
update_nlist) |
| 414 |
#endif |
| 415 |
endif |
| 416 |
|
| 417 |
if (update_nlist) then |
| 418 |
!! save current configuration and construct neighbor list |
| 419 |
#ifdef IS_MPI |
| 420 |
call saveNeighborList(nGroupsInRow, q_group_row) |
| 421 |
#else |
| 422 |
call saveNeighborList(nGroups, q_group) |
| 423 |
#endif |
| 424 |
neighborListSize = size(list) |
| 425 |
nlist = 0 |
| 426 |
endif |
| 427 |
|
| 428 |
istart = 1 |
| 429 |
#ifdef IS_MPI |
| 430 |
iend = nGroupsInRow |
| 431 |
#else |
| 432 |
iend = nGroups - 1 |
| 433 |
#endif |
| 434 |
outer: do i = istart, iend |
| 435 |
|
| 436 |
if (update_nlist) point(i) = nlist + 1 |
| 437 |
|
| 438 |
n_in_i = groupStartRow(i+1) - groupStartRow(i) |
| 439 |
|
| 440 |
if (update_nlist) then |
| 441 |
#ifdef IS_MPI |
| 442 |
jstart = 1 |
| 443 |
jend = nGroupsInCol |
| 444 |
#else |
| 445 |
jstart = i+1 |
| 446 |
jend = nGroups |
| 447 |
#endif |
| 448 |
else |
| 449 |
jstart = point(i) |
| 450 |
jend = point(i+1) - 1 |
| 451 |
! make sure group i has neighbors |
| 452 |
if (jstart .gt. jend) cycle outer |
| 453 |
endif |
| 454 |
|
| 455 |
do jnab = jstart, jend |
| 456 |
if (update_nlist) then |
| 457 |
j = jnab |
| 458 |
else |
| 459 |
j = list(jnab) |
| 460 |
endif |
| 461 |
|
| 462 |
#ifdef IS_MPI |
| 463 |
me_j = atid_col(j) |
| 464 |
call get_interatomic_vector(q_group_Row(:,i), & |
| 465 |
q_group_Col(:,j), d_grp, rgrpsq) |
| 466 |
#else |
| 467 |
me_j = atid(j) |
| 468 |
call get_interatomic_vector(q_group(:,i), & |
| 469 |
q_group(:,j), d_grp, rgrpsq) |
| 470 |
#endif |
| 471 |
|
| 472 |
if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rListsq) then |
| 473 |
if (update_nlist) then |
| 474 |
nlist = nlist + 1 |
| 475 |
|
| 476 |
if (nlist > neighborListSize) then |
| 477 |
#ifdef IS_MPI |
| 478 |
call expandNeighborList(nGroupsInRow, listerror) |
| 479 |
#else |
| 480 |
call expandNeighborList(nGroups, listerror) |
| 481 |
#endif |
| 482 |
if (listerror /= 0) then |
| 483 |
error = -1 |
| 484 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
| 485 |
return |
| 486 |
end if |
| 487 |
neighborListSize = size(list) |
| 488 |
endif |
| 489 |
|
| 490 |
list(nlist) = j |
| 491 |
endif |
| 492 |
|
| 493 |
if (rgrpsq < gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCutsq) then |
| 494 |
|
| 495 |
rCut = gtypeCutoffMap(groupToGtypeRow(i),groupToGtypeCol(j))%rCut |
| 496 |
if (loop .eq. PAIR_LOOP) then |
| 497 |
vij = 0.0_dp |
| 498 |
fij(1) = 0.0_dp |
| 499 |
fij(2) = 0.0_dp |
| 500 |
fij(3) = 0.0_dp |
| 501 |
endif |
| 502 |
|
| 503 |
call get_switch(rgrpsq, sw, dswdr,rgrp, in_switching_region) |
| 504 |
|
| 505 |
n_in_j = groupStartCol(j+1) - groupStartCol(j) |
| 506 |
|
| 507 |
do ia = groupStartRow(i), groupStartRow(i+1)-1 |
| 508 |
|
| 509 |
atom1 = groupListRow(ia) |
| 510 |
|
| 511 |
inner: do jb = groupStartCol(j), groupStartCol(j+1)-1 |
| 512 |
|
| 513 |
atom2 = groupListCol(jb) |
| 514 |
|
| 515 |
if (skipThisPair(atom1, atom2)) cycle inner |
| 516 |
|
| 517 |
if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then |
| 518 |
d_atm(1) = d_grp(1) |
| 519 |
d_atm(2) = d_grp(2) |
| 520 |
d_atm(3) = d_grp(3) |
| 521 |
ratmsq = rgrpsq |
| 522 |
else |
| 523 |
#ifdef IS_MPI |
| 524 |
call get_interatomic_vector(q_Row(:,atom1), & |
| 525 |
q_Col(:,atom2), d_atm, ratmsq) |
| 526 |
#else |
| 527 |
call get_interatomic_vector(q(:,atom1), & |
| 528 |
q(:,atom2), d_atm, ratmsq) |
| 529 |
#endif |
| 530 |
endif |
| 531 |
|
| 532 |
topoDist = getTopoDistance(atom1, atom2) |
| 533 |
|
| 534 |
if (loop .eq. PREPAIR_LOOP) then |
| 535 |
#ifdef IS_MPI |
| 536 |
call do_prepair(atom1, atom2, ratmsq, d_atm, sw, & |
| 537 |
rgrpsq, d_grp, rCut, & |
| 538 |
eFrame, A, f, t, pot_local) |
| 539 |
#else |
| 540 |
call do_prepair(atom1, atom2, ratmsq, d_atm, sw, & |
| 541 |
rgrpsq, d_grp, rCut, & |
| 542 |
eFrame, A, f, t, pot) |
| 543 |
#endif |
| 544 |
else |
| 545 |
#ifdef IS_MPI |
| 546 |
call f_do_pair(atom1, atom2, ratmsq, d_atm, sw, & |
| 547 |
eFrame, A, f, t, pot_local, particle_pot, vpair, & |
| 548 |
fpair, d_grp, rgrp, rCut, topoDist) |
| 549 |
! particle_pot will be accumulated from row & column |
| 550 |
! arrays later |
| 551 |
#else |
| 552 |
call f_do_pair(atom1, atom2, ratmsq, d_atm, sw, & |
| 553 |
eFrame, A, f, t, pot, particle_pot, vpair, & |
| 554 |
fpair, d_grp, rgrp, rCut, topoDist) |
| 555 |
#endif |
| 556 |
vij = vij + vpair |
| 557 |
fij(1) = fij(1) + fpair(1) |
| 558 |
fij(2) = fij(2) + fpair(2) |
| 559 |
fij(3) = fij(3) + fpair(3) |
| 560 |
call add_stress_tensor(d_atm, fpair, tau) |
| 561 |
endif |
| 562 |
enddo inner |
| 563 |
enddo |
| 564 |
|
| 565 |
if (loop .eq. PAIR_LOOP) then |
| 566 |
if (in_switching_region) then |
| 567 |
swderiv = vij*dswdr/rgrp |
| 568 |
fg = swderiv*d_grp |
| 569 |
|
| 570 |
fij(1) = fij(1) + fg(1) |
| 571 |
fij(2) = fij(2) + fg(2) |
| 572 |
fij(3) = fij(3) + fg(3) |
| 573 |
|
| 574 |
if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then |
| 575 |
call add_stress_tensor(d_atm, fg, tau) |
| 576 |
endif |
| 577 |
|
| 578 |
do ia=groupStartRow(i), groupStartRow(i+1)-1 |
| 579 |
atom1=groupListRow(ia) |
| 580 |
mf = mfactRow(atom1) |
| 581 |
! fg is the force on atom ia due to cutoff group's |
| 582 |
! presence in switching region |
| 583 |
fg = swderiv*d_grp*mf |
| 584 |
#ifdef IS_MPI |
| 585 |
f_Row(1,atom1) = f_Row(1,atom1) + fg(1) |
| 586 |
f_Row(2,atom1) = f_Row(2,atom1) + fg(2) |
| 587 |
f_Row(3,atom1) = f_Row(3,atom1) + fg(3) |
| 588 |
#else |
| 589 |
f(1,atom1) = f(1,atom1) + fg(1) |
| 590 |
f(2,atom1) = f(2,atom1) + fg(2) |
| 591 |
f(3,atom1) = f(3,atom1) + fg(3) |
| 592 |
#endif |
| 593 |
if (n_in_i .gt. 1) then |
| 594 |
if (SIM_uses_AtomicVirial) then |
| 595 |
! find the distance between the atom |
| 596 |
! and the center of the cutoff group: |
| 597 |
#ifdef IS_MPI |
| 598 |
call get_interatomic_vector(q_Row(:,atom1), & |
| 599 |
q_group_Row(:,i), dag, rag) |
| 600 |
#else |
| 601 |
call get_interatomic_vector(q(:,atom1), & |
| 602 |
q_group(:,i), dag, rag) |
| 603 |
#endif |
| 604 |
call add_stress_tensor(dag,fg,tau) |
| 605 |
endif |
| 606 |
endif |
| 607 |
enddo |
| 608 |
|
| 609 |
do jb=groupStartCol(j), groupStartCol(j+1)-1 |
| 610 |
atom2=groupListCol(jb) |
| 611 |
mf = mfactCol(atom2) |
| 612 |
! fg is the force on atom jb due to cutoff group's |
| 613 |
! presence in switching region |
| 614 |
fg = -swderiv*d_grp*mf |
| 615 |
#ifdef IS_MPI |
| 616 |
f_Col(1,atom2) = f_Col(1,atom2) + fg(1) |
| 617 |
f_Col(2,atom2) = f_Col(2,atom2) + fg(2) |
| 618 |
f_Col(3,atom2) = f_Col(3,atom2) + fg(3) |
| 619 |
#else |
| 620 |
f(1,atom2) = f(1,atom2) + fg(1) |
| 621 |
f(2,atom2) = f(2,atom2) + fg(2) |
| 622 |
f(3,atom2) = f(3,atom2) + fg(3) |
| 623 |
#endif |
| 624 |
if (n_in_j .gt. 1) then |
| 625 |
if (SIM_uses_AtomicVirial) then |
| 626 |
! find the distance between the atom |
| 627 |
! and the center of the cutoff group: |
| 628 |
#ifdef IS_MPI |
| 629 |
call get_interatomic_vector(q_Col(:,atom2), & |
| 630 |
q_group_Col(:,j), dag, rag) |
| 631 |
#else |
| 632 |
call get_interatomic_vector(q(:,atom2), & |
| 633 |
q_group(:,j), dag, rag) |
| 634 |
#endif |
| 635 |
call add_stress_tensor(dag,fg,tau) |
| 636 |
endif |
| 637 |
endif |
| 638 |
enddo |
| 639 |
endif |
| 640 |
!if (.not.SIM_uses_AtomicVirial) then |
| 641 |
! call add_stress_tensor(d_grp, fij, tau) |
| 642 |
!endif |
| 643 |
endif |
| 644 |
endif |
| 645 |
endif |
| 646 |
enddo |
| 647 |
|
| 648 |
enddo outer |
| 649 |
|
| 650 |
if (update_nlist) then |
| 651 |
#ifdef IS_MPI |
| 652 |
point(nGroupsInRow + 1) = nlist + 1 |
| 653 |
#else |
| 654 |
point(nGroups) = nlist + 1 |
| 655 |
#endif |
| 656 |
if (loop .eq. PREPAIR_LOOP) then |
| 657 |
! we just did the neighbor list update on the first |
| 658 |
! pass, so we don't need to do it |
| 659 |
! again on the second pass |
| 660 |
update_nlist = .false. |
| 661 |
endif |
| 662 |
endif |
| 663 |
|
| 664 |
if (loop .eq. PREPAIR_LOOP) then |
| 665 |
#ifdef IS_MPI |
| 666 |
call do_preforce(nlocal, pot_local, particle_pot) |
| 667 |
#else |
| 668 |
call do_preforce(nlocal, pot, particle_pot) |
| 669 |
#endif |
| 670 |
endif |
| 671 |
|
| 672 |
enddo |
| 673 |
|
| 674 |
!! Do timing |
| 675 |
#ifdef PROFILE |
| 676 |
call cpu_time(forceTimeFinal) |
| 677 |
forceTime = forceTime + forceTimeFinal - forceTimeInitial |
| 678 |
#endif |
| 679 |
|
| 680 |
#ifdef IS_MPI |
| 681 |
!!distribute forces |
| 682 |
|
| 683 |
f_temp = 0.0_dp |
| 684 |
call scatter(f_Row,f_temp,plan_atom_row_3d) |
| 685 |
do i = 1,nlocal |
| 686 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
| 687 |
end do |
| 688 |
|
| 689 |
f_temp = 0.0_dp |
| 690 |
call scatter(f_Col,f_temp,plan_atom_col_3d) |
| 691 |
do i = 1,nlocal |
| 692 |
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
| 693 |
end do |
| 694 |
|
| 695 |
if (FF_UsesDirectionalAtoms() .and. SIM_uses_DirectionalAtoms) then |
| 696 |
t_temp = 0.0_dp |
| 697 |
call scatter(t_Row,t_temp,plan_atom_row_3d) |
| 698 |
do i = 1,nlocal |
| 699 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
| 700 |
end do |
| 701 |
t_temp = 0.0_dp |
| 702 |
call scatter(t_Col,t_temp,plan_atom_col_3d) |
| 703 |
|
| 704 |
do i = 1,nlocal |
| 705 |
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
| 706 |
end do |
| 707 |
endif |
| 708 |
|
| 709 |
! scatter/gather pot_row into the members of my column |
| 710 |
do i = 1,LR_POT_TYPES |
| 711 |
call scatter(pot_Row(i,:), pot_Temp(i,:), plan_atom_row) |
| 712 |
end do |
| 713 |
! scatter/gather pot_local into all other procs |
| 714 |
! add resultant to get total pot |
| 715 |
do i = 1, nlocal |
| 716 |
pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES) & |
| 717 |
+ pot_Temp(1:LR_POT_TYPES,i) |
| 718 |
enddo |
| 719 |
|
| 720 |
do i = 1,LR_POT_TYPES |
| 721 |
particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal) |
| 722 |
enddo |
| 723 |
|
| 724 |
pot_Temp = 0.0_DP |
| 725 |
|
| 726 |
do i = 1,LR_POT_TYPES |
| 727 |
call scatter(pot_Col(i,:), pot_Temp(i,:), plan_atom_col) |
| 728 |
end do |
| 729 |
|
| 730 |
do i = 1, nlocal |
| 731 |
pot_local(1:LR_POT_TYPES) = pot_local(1:LR_POT_TYPES)& |
| 732 |
+ pot_Temp(1:LR_POT_TYPES,i) |
| 733 |
enddo |
| 734 |
|
| 735 |
do i = 1,LR_POT_TYPES |
| 736 |
particle_pot(1:nlocal) = particle_pot(1:nlocal) + pot_Temp(i,1:nlocal) |
| 737 |
enddo |
| 738 |
|
| 739 |
ppot_Temp = 0.0_DP |
| 740 |
|
| 741 |
call scatter(ppot_Row(:), ppot_Temp(:), plan_atom_row) |
| 742 |
do i = 1, nlocal |
| 743 |
particle_pot(i) = particle_pot(i) + ppot_Temp(i) |
| 744 |
enddo |
| 745 |
|
| 746 |
ppot_Temp = 0.0_DP |
| 747 |
|
| 748 |
call scatter(ppot_Col(:), ppot_Temp(:), plan_atom_col) |
| 749 |
do i = 1, nlocal |
| 750 |
particle_pot(i) = particle_pot(i) + ppot_Temp(i) |
| 751 |
enddo |
| 752 |
|
| 753 |
#endif |
| 754 |
|
| 755 |
if (SIM_requires_postpair_calc) then |
| 756 |
do i = 1, nlocal |
| 757 |
|
| 758 |
do i1 = 1, nSkipsForLocalAtom(i) |
| 759 |
j = skipsForLocalAtom(i, i1) |
| 760 |
|
| 761 |
! prevent overcounting the skips |
| 762 |
if (i.lt.j) then |
| 763 |
|
| 764 |
call get_interatomic_vector(q(:,i), q(:,j), d_atm, ratmsq) |
| 765 |
rVal = sqrt(ratmsq) |
| 766 |
call get_switch(ratmsq, sw, dswdr, rVal,in_switching_region) |
| 767 |
#ifdef IS_MPI |
| 768 |
call do_skip_correction(c_idents_local(i), c_idents_local(j), & |
| 769 |
d_atm, rVal, skipped_charge(i), skipped_charge(j), sw, & |
| 770 |
1.0_dp, pot_local(ELECTROSTATIC_POT), vpair, f, t(:,i), t(:,j)) |
| 771 |
# else |
| 772 |
call do_skip_correction(c_idents_local(i), c_idents_local(j), & |
| 773 |
d_atm, rVal, skipped_charge(i), skipped_charge(j), sw, & |
| 774 |
1.0_dp, pot(ELECTROSTATIC_POT), vpair, f, t(:,i), t(:,j)) |
| 775 |
#endif |
| 776 |
endif |
| 777 |
enddo |
| 778 |
enddo |
| 779 |
|
| 780 |
do i = 1, nlocal |
| 781 |
|
| 782 |
#ifdef IS_MPI |
| 783 |
call do_self_correction(c_idents_local(i), eFrame(:,i), & |
| 784 |
skipped_charge(i), pot_local(ELECTROSTATIC_POT), t(:,i)) |
| 785 |
#else |
| 786 |
call do_self_correction(c_idents_local(i), eFrame(:,i), & |
| 787 |
skipped_charge(i), pot(ELECTROSTATIC_POT), t(:,i)) |
| 788 |
#endif |
| 789 |
enddo |
| 790 |
endif |
| 791 |
|
| 792 |
#ifdef IS_MPI |
| 793 |
#ifdef SINGLE_PRECISION |
| 794 |
call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_real,mpi_sum, & |
| 795 |
mpi_comm_world,mpi_err) |
| 796 |
#else |
| 797 |
call mpi_allreduce(pot_local, pot, LR_POT_TYPES,mpi_double_precision, & |
| 798 |
mpi_sum, mpi_comm_world,mpi_err) |
| 799 |
#endif |
| 800 |
#endif |
| 801 |
|
| 802 |
end subroutine do_force_loop |
| 803 |
|
| 804 |
subroutine f_do_pair(i, j, rijsq, d, sw, & |
| 805 |
eFrame, A, f, t, pot, particle_pot, vpair, & |
| 806 |
fpair, d_grp, r_grp, rCut, topoDist) |
| 807 |
|
| 808 |
real( kind = dp ) :: vpair, sw |
| 809 |
real( kind = dp ), dimension(LR_POT_TYPES) :: pot, pairpot |
| 810 |
real( kind = dp ), dimension(nLocal) :: particle_pot |
| 811 |
real( kind = dp ), dimension(3) :: fpair |
| 812 |
real( kind = dp ), dimension(nLocal) :: mfact |
| 813 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
| 814 |
real( kind = dp ), dimension(9,nLocal) :: A |
| 815 |
real( kind = dp ), dimension(3,nLocal) :: f |
| 816 |
real( kind = dp ), dimension(3,nLocal) :: t |
| 817 |
|
| 818 |
integer, intent(in) :: i, j |
| 819 |
real ( kind = dp ), intent(inout) :: rijsq |
| 820 |
real ( kind = dp ), intent(inout) :: r_grp |
| 821 |
real ( kind = dp ), intent(inout) :: d(3) |
| 822 |
real ( kind = dp ), intent(inout) :: d_grp(3) |
| 823 |
real ( kind = dp ), intent(inout) :: rCut |
| 824 |
integer, intent(inout) :: topoDist |
| 825 |
real ( kind = dp ) :: r, pair_pot, vdwMult, electroMult |
| 826 |
real ( kind = dp ) :: a_k, b_k, c_k, d_k, dx |
| 827 |
|
| 828 |
real( kind = dp), dimension(3) :: f1, t1, t2 |
| 829 |
real( kind = dp), dimension(9) :: A1, A2, eF1, eF2 |
| 830 |
real( kind = dp) :: dfrhodrho_i, dfrhodrho_j |
| 831 |
real( kind = dp) :: rho_i, rho_j |
| 832 |
real( kind = dp) :: fshift_i, fshift_j |
| 833 |
real( kind = dp) :: p_vdw, p_elect, p_hb, p_met |
| 834 |
integer :: id1, id2, idx |
| 835 |
integer :: k |
| 836 |
integer :: c_ident_i, c_ident_j |
| 837 |
|
| 838 |
integer :: iHash |
| 839 |
|
| 840 |
r = sqrt(rijsq) |
| 841 |
|
| 842 |
vpair = 0.0_dp |
| 843 |
fpair(1:3) = 0.0_dp |
| 844 |
|
| 845 |
p_vdw = 0.0 |
| 846 |
p_elect = 0.0 |
| 847 |
p_hb = 0.0 |
| 848 |
p_met = 0.0 |
| 849 |
|
| 850 |
f1(1:3) = 0.0 |
| 851 |
t1(1:3) = 0.0 |
| 852 |
t2(1:3) = 0.0 |
| 853 |
|
| 854 |
#ifdef IS_MPI |
| 855 |
c_ident_i = c_idents_row(i) |
| 856 |
c_ident_j = c_idents_col(j) |
| 857 |
|
| 858 |
A1(:) = A_Row(:,i) |
| 859 |
A2(:) = A_Col(:,j) |
| 860 |
eF1(:) = eFrame_Row(:,i) |
| 861 |
eF2(:) = eFrame_Col(:,j) |
| 862 |
|
| 863 |
dfrhodrho_i = dfrhodrho_row(i) |
| 864 |
dfrhodrho_j = dfrhodrho_col(j) |
| 865 |
rho_i = rho_row(i) |
| 866 |
rho_j = rho_col(j) |
| 867 |
#else |
| 868 |
c_ident_i = c_idents_local(i) |
| 869 |
c_ident_j = c_idents_local(j) |
| 870 |
|
| 871 |
A1(:) = A(:,i) |
| 872 |
A2(:) = A(:,j) |
| 873 |
eF1(:) = eFrame(:,i) |
| 874 |
eF2(:) = eFrame(:,j) |
| 875 |
|
| 876 |
dfrhodrho_i = dfrhodrho(i) |
| 877 |
dfrhodrho_j = dfrhodrho(j) |
| 878 |
rho_i = rho(i) |
| 879 |
rho_j = rho(j) |
| 880 |
#endif |
| 881 |
|
| 882 |
vdwMult = vdwScale(topoDist) |
| 883 |
electroMult = electrostaticScale(topoDist) |
| 884 |
|
| 885 |
call doPairInteraction(c_ident_i, c_ident_j, d, r, rijsq, sw, vpair, & |
| 886 |
vdwMult, electroMult, A1, A2, eF1, eF2, & |
| 887 |
pairpot, f1, t1, t2, & |
| 888 |
rho_i, rho_j, dfrhodrho_i, dfrhodrho_j, fshift_i, fshift_j) |
| 889 |
|
| 890 |
#ifdef IS_MPI |
| 891 |
id1 = AtomRowToGlobal(i) |
| 892 |
id2 = AtomColToGlobal(j) |
| 893 |
|
| 894 |
pot_row(VDW_POT,i) = pot_row(VDW_POT,i) + 0.5*pairpot(VDW_POT) |
| 895 |
pot_col(VDW_POT,j) = pot_col(VDW_POT,j) + 0.5*pairpot(VDW_POT) |
| 896 |
pot_row(ELECTROSTATIC_POT,i) = pot_row(ELECTROSTATIC_POT,i) + 0.5*pairpot(ELECTROSTATIC_POT) |
| 897 |
pot_col(ELECTROSTATIC_POT,j) = pot_col(ELECTROSTATIC_POT,j) + 0.5*pairpot(ELECTROSTATIC_POT) |
| 898 |
pot_row(HB_POT,i) = pot_row(HB_POT,i) + 0.5*pairpot(HB_POT) |
| 899 |
pot_col(HB_POT,j) = pot_col(HB_POT,j) + 0.5*pairpot(HB_POT) |
| 900 |
pot_Row(METALLIC_POT,i) = pot_Row(METALLIC_POT,i) + 0.5*pairpot(METALLIC_POT) |
| 901 |
pot_Col(METALLIC_POT,j) = pot_Col(METALLIC_POT,j) + 0.5*pairpot(METALLIC_POT) |
| 902 |
|
| 903 |
do idx = 1, 3 |
| 904 |
f_Row(idx,i) = f_Row(idx,i) + f1(idx) |
| 905 |
f_Col(idx,j) = f_Col(idx,j) - f1(idx) |
| 906 |
|
| 907 |
t_Row(idx,i) = t_Row(idx,i) + t1(idx) |
| 908 |
t_Col(idx,j) = t_Col(idx,j) + t2(idx) |
| 909 |
enddo |
| 910 |
! particle_pot is the difference between the full potential |
| 911 |
! and the full potential without the presence of a particular |
| 912 |
! particle (atom1). |
| 913 |
! |
| 914 |
! This reduces the density at other particle locations, so |
| 915 |
! we need to recompute the density at atom2 assuming atom1 |
| 916 |
! didn't contribute. This then requires recomputing the |
| 917 |
! density functional for atom2 as well. |
| 918 |
! |
| 919 |
! Most of the particle_pot heavy lifting comes from the |
| 920 |
! pair interaction, and will be handled by vpair. Parallel version. |
| 921 |
|
| 922 |
if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then |
| 923 |
ppot_row(i) = ppot_row(i) - frho_row(j) + fshift_j |
| 924 |
ppot_col(j) = ppot_col(j) - frho_col(i) + fshift_i |
| 925 |
end if |
| 926 |
|
| 927 |
#else |
| 928 |
id1 = i |
| 929 |
id2 = j |
| 930 |
|
| 931 |
pot(VDW_POT) = pot(VDW_POT) + pairpot(VDW_POT) |
| 932 |
pot(ELECTROSTATIC_POT) = pot(ELECTROSTATIC_POT) + pairpot(ELECTROSTATIC_POT) |
| 933 |
pot(HB_POT) = pot(HB_POT) + pairpot(HB_POT) |
| 934 |
pot(METALLIC_POT) = pot(METALLIC_POT) + pairpot(METALLIC_POT) |
| 935 |
|
| 936 |
do idx = 1, 3 |
| 937 |
f(idx,i) = f(idx,i) + f1(idx) |
| 938 |
f(idx,j) = f(idx,j) - f1(idx) |
| 939 |
|
| 940 |
t(idx,i) = t(idx,i) + t1(idx) |
| 941 |
t(idx,j) = t(idx,j) + t2(idx) |
| 942 |
enddo |
| 943 |
! particle_pot is the difference between the full potential |
| 944 |
! and the full potential without the presence of a particular |
| 945 |
! particle (atom1). |
| 946 |
! |
| 947 |
! This reduces the density at other particle locations, so |
| 948 |
! we need to recompute the density at atom2 assuming atom1 |
| 949 |
! didn't contribute. This then requires recomputing the |
| 950 |
! density functional for atom2 as well. |
| 951 |
! |
| 952 |
! Most of the particle_pot heavy lifting comes from the |
| 953 |
! pair interaction, and will be handled by vpair. NonParallel version. |
| 954 |
|
| 955 |
if ( (iand(iHash, EAM_PAIR).ne.0) .or. (iand(iHash, SC_PAIR).ne.0) ) then |
| 956 |
particle_pot(i) = particle_pot(i) - frho(j) + fshift_j |
| 957 |
particle_pot(j) = particle_pot(j) - frho(i) + fshift_i |
| 958 |
end if |
| 959 |
|
| 960 |
|
| 961 |
#endif |
| 962 |
|
| 963 |
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
| 964 |
|
| 965 |
fpair(1) = fpair(1) + f1(1) |
| 966 |
fpair(2) = fpair(2) + f1(2) |
| 967 |
fpair(3) = fpair(3) + f1(3) |
| 968 |
|
| 969 |
endif |
| 970 |
end subroutine f_do_pair |
| 971 |
|
| 972 |
subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, rCut, & |
| 973 |
eFrame, A, f, t, pot) |
| 974 |
|
| 975 |
real( kind = dp ) :: sw |
| 976 |
real( kind = dp ), dimension(LR_POT_TYPES) :: pot |
| 977 |
real( kind = dp ), dimension(9,nLocal) :: eFrame |
| 978 |
real (kind=dp), dimension(9,nLocal) :: A |
| 979 |
real (kind=dp), dimension(3,nLocal) :: f |
| 980 |
real (kind=dp), dimension(3,nLocal) :: t |
| 981 |
|
| 982 |
integer, intent(in) :: i, j |
| 983 |
real ( kind = dp ), intent(inout) :: rijsq, rcijsq, rCut |
| 984 |
real ( kind = dp ) :: r, rc |
| 985 |
real ( kind = dp ), intent(inout) :: d(3), dc(3) |
| 986 |
real ( kind = dp ) :: rho_i_at_j, rho_j_at_i |
| 987 |
integer :: c_ident_i, c_ident_j |
| 988 |
|
| 989 |
r = sqrt(rijsq) |
| 990 |
|
| 991 |
#ifdef IS_MPI |
| 992 |
c_ident_i = c_idents_row(i) |
| 993 |
c_ident_j = c_idents_col(j) |
| 994 |
#else |
| 995 |
c_ident_i = c_idents_local(i) |
| 996 |
c_ident_j = c_idents_local(j) |
| 997 |
#endif |
| 998 |
rho_i_at_j = 0.0_dp |
| 999 |
rho_j_at_i = 0.0_dp |
| 1000 |
|
| 1001 |
call doPrepairInteraction(c_ident_i, c_ident_j, r, & |
| 1002 |
rho_i_at_j, rho_j_at_i) |
| 1003 |
|
| 1004 |
#ifdef IS_MPI |
| 1005 |
rho_col(j) = rho_col(j) + rho_i_at_j |
| 1006 |
rho_row(i) = rho_row(i) + rho_j_at_i |
| 1007 |
#else |
| 1008 |
rho(j) = rho(j) + rho_i_at_j |
| 1009 |
rho(i) = rho(i) + rho_j_at_i |
| 1010 |
#endif |
| 1011 |
|
| 1012 |
end subroutine do_prepair |
| 1013 |
|
| 1014 |
|
| 1015 |
subroutine do_preforce(nlocal, pot, particle_pot) |
| 1016 |
integer :: nlocal |
| 1017 |
real( kind = dp ),dimension(LR_POT_TYPES) :: pot |
| 1018 |
real( kind = dp ),dimension(nlocal) :: particle_pot |
| 1019 |
integer :: sc_err = 0 |
| 1020 |
integer :: atom, c_ident1 |
| 1021 |
|
| 1022 |
if ((FF_uses_EAM .and. SIM_uses_EAM) .or. (FF_uses_SC .and. SIM_uses_SC)) then |
| 1023 |
|
| 1024 |
#ifdef IS_MPI |
| 1025 |
call scatter(rho_row,rho,plan_atom_row,sc_err) |
| 1026 |
if (sc_err /= 0 ) then |
| 1027 |
call handleError("do_preforce()", "Error scattering rho_row into rho") |
| 1028 |
endif |
| 1029 |
call scatter(rho_col,rho_tmp,plan_atom_col,sc_err) |
| 1030 |
if (sc_err /= 0 ) then |
| 1031 |
call handleError("do_preforce()", "Error scattering rho_col into rho") |
| 1032 |
endif |
| 1033 |
rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal) |
| 1034 |
#endif |
| 1035 |
|
| 1036 |
|
| 1037 |
do atom = 1, nlocal |
| 1038 |
c_ident1 = c_idents_local(atom) |
| 1039 |
|
| 1040 |
call doPreforceInteraction(c_ident1, rho(atom), frho(atom), dfrhodrho(atom)) |
| 1041 |
pot(METALLIC_POT) = pot(METALLIC_POT) + frho(atom) |
| 1042 |
particle_pot(atom) = particle_pot(atom) + frho(atom) |
| 1043 |
end do |
| 1044 |
|
| 1045 |
#ifdef IS_MPI |
| 1046 |
!! communicate f(rho) and derivatives back into row and column arrays |
| 1047 |
call gather(frho,frho_row,plan_atom_row, sc_err) |
| 1048 |
if (sc_err /= 0) then |
| 1049 |
call handleError("do_preforce()","MPI gather frho_row failure") |
| 1050 |
endif |
| 1051 |
call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err) |
| 1052 |
if (sc_err /= 0) then |
| 1053 |
call handleError("do_preforce()","MPI gather dfrhodrho_row failure") |
| 1054 |
endif |
| 1055 |
call gather(frho,frho_col,plan_atom_col, sc_err) |
| 1056 |
if (sc_err /= 0) then |
| 1057 |
call handleError("do_preforce()","MPI gather frho_col failure") |
| 1058 |
endif |
| 1059 |
call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err) |
| 1060 |
if (sc_err /= 0) then |
| 1061 |
call handleError("do_preforce()","MPI gather dfrhodrho_col failure") |
| 1062 |
endif |
| 1063 |
#endif |
| 1064 |
|
| 1065 |
end if |
| 1066 |
end subroutine do_preforce |
| 1067 |
|
| 1068 |
|
| 1069 |
subroutine get_interatomic_vector(q_i, q_j, d, r_sq) |
| 1070 |
|
| 1071 |
real (kind = dp), dimension(3) :: q_i |
| 1072 |
real (kind = dp), dimension(3) :: q_j |
| 1073 |
real ( kind = dp ), intent(out) :: r_sq |
| 1074 |
real( kind = dp ) :: d(3), scaled(3) |
| 1075 |
integer i |
| 1076 |
|
| 1077 |
d(1) = q_j(1) - q_i(1) |
| 1078 |
d(2) = q_j(2) - q_i(2) |
| 1079 |
d(3) = q_j(3) - q_i(3) |
| 1080 |
|
| 1081 |
! Wrap back into periodic box if necessary |
| 1082 |
if ( SIM_uses_PBC ) then |
| 1083 |
|
| 1084 |
if( .not.boxIsOrthorhombic ) then |
| 1085 |
! calc the scaled coordinates. |
| 1086 |
! scaled = matmul(HmatInv, d) |
| 1087 |
|
| 1088 |
scaled(1) = HmatInv(1,1)*d(1) + HmatInv(1,2)*d(2) + HmatInv(1,3)*d(3) |
| 1089 |
scaled(2) = HmatInv(2,1)*d(1) + HmatInv(2,2)*d(2) + HmatInv(2,3)*d(3) |
| 1090 |
scaled(3) = HmatInv(3,1)*d(1) + HmatInv(3,2)*d(2) + HmatInv(3,3)*d(3) |
| 1091 |
|
| 1092 |
! wrap the scaled coordinates |
| 1093 |
|
| 1094 |
scaled(1) = scaled(1) - anint(scaled(1), kind=dp) |
| 1095 |
scaled(2) = scaled(2) - anint(scaled(2), kind=dp) |
| 1096 |
scaled(3) = scaled(3) - anint(scaled(3), kind=dp) |
| 1097 |
|
| 1098 |
! calc the wrapped real coordinates from the wrapped scaled |
| 1099 |
! coordinates |
| 1100 |
! d = matmul(Hmat,scaled) |
| 1101 |
d(1)= Hmat(1,1)*scaled(1) + Hmat(1,2)*scaled(2) + Hmat(1,3)*scaled(3) |
| 1102 |
d(2)= Hmat(2,1)*scaled(1) + Hmat(2,2)*scaled(2) + Hmat(2,3)*scaled(3) |
| 1103 |
d(3)= Hmat(3,1)*scaled(1) + Hmat(3,2)*scaled(2) + Hmat(3,3)*scaled(3) |
| 1104 |
|
| 1105 |
else |
| 1106 |
! calc the scaled coordinates. |
| 1107 |
|
| 1108 |
scaled(1) = d(1) * HmatInv(1,1) |
| 1109 |
scaled(2) = d(2) * HmatInv(2,2) |
| 1110 |
scaled(3) = d(3) * HmatInv(3,3) |
| 1111 |
|
| 1112 |
! wrap the scaled coordinates |
| 1113 |
|
| 1114 |
scaled(1) = scaled(1) - anint(scaled(1), kind=dp) |
| 1115 |
scaled(2) = scaled(2) - anint(scaled(2), kind=dp) |
| 1116 |
scaled(3) = scaled(3) - anint(scaled(3), kind=dp) |
| 1117 |
|
| 1118 |
! calc the wrapped real coordinates from the wrapped scaled |
| 1119 |
! coordinates |
| 1120 |
|
| 1121 |
d(1) = scaled(1)*Hmat(1,1) |
| 1122 |
d(2) = scaled(2)*Hmat(2,2) |
| 1123 |
d(3) = scaled(3)*Hmat(3,3) |
| 1124 |
|
| 1125 |
endif |
| 1126 |
|
| 1127 |
endif |
| 1128 |
|
| 1129 |
r_sq = d(1)*d(1) + d(2)*d(2) + d(3)*d(3) |
| 1130 |
|
| 1131 |
end subroutine get_interatomic_vector |
| 1132 |
|
| 1133 |
subroutine zero_work_arrays() |
| 1134 |
|
| 1135 |
#ifdef IS_MPI |
| 1136 |
|
| 1137 |
q_Row = 0.0_dp |
| 1138 |
q_Col = 0.0_dp |
| 1139 |
|
| 1140 |
q_group_Row = 0.0_dp |
| 1141 |
q_group_Col = 0.0_dp |
| 1142 |
|
| 1143 |
eFrame_Row = 0.0_dp |
| 1144 |
eFrame_Col = 0.0_dp |
| 1145 |
|
| 1146 |
A_Row = 0.0_dp |
| 1147 |
A_Col = 0.0_dp |
| 1148 |
|
| 1149 |
f_Row = 0.0_dp |
| 1150 |
f_Col = 0.0_dp |
| 1151 |
f_Temp = 0.0_dp |
| 1152 |
|
| 1153 |
t_Row = 0.0_dp |
| 1154 |
t_Col = 0.0_dp |
| 1155 |
t_Temp = 0.0_dp |
| 1156 |
|
| 1157 |
pot_Row = 0.0_dp |
| 1158 |
pot_Col = 0.0_dp |
| 1159 |
pot_Temp = 0.0_dp |
| 1160 |
ppot_Temp = 0.0_dp |
| 1161 |
|
| 1162 |
frho_row = 0.0_dp |
| 1163 |
frho_col = 0.0_dp |
| 1164 |
rho_row = 0.0_dp |
| 1165 |
rho_col = 0.0_dp |
| 1166 |
rho_tmp = 0.0_dp |
| 1167 |
dfrhodrho_row = 0.0_dp |
| 1168 |
dfrhodrho_col = 0.0_dp |
| 1169 |
|
| 1170 |
#endif |
| 1171 |
rho = 0.0_dp |
| 1172 |
frho = 0.0_dp |
| 1173 |
dfrhodrho = 0.0_dp |
| 1174 |
|
| 1175 |
end subroutine zero_work_arrays |
| 1176 |
|
| 1177 |
function skipThisPair(atom1, atom2) result(skip_it) |
| 1178 |
integer, intent(in) :: atom1 |
| 1179 |
integer, intent(in), optional :: atom2 |
| 1180 |
logical :: skip_it |
| 1181 |
integer :: unique_id_1, unique_id_2 |
| 1182 |
integer :: i |
| 1183 |
|
| 1184 |
skip_it = .false. |
| 1185 |
|
| 1186 |
!! there are a number of reasons to skip a pair or a particle |
| 1187 |
!! mostly we do this to exclude atoms who are involved in short |
| 1188 |
!! range interactions (bonds, bends, torsions), but we also need |
| 1189 |
!! to exclude some overcounted interactions that result from |
| 1190 |
!! the parallel decomposition |
| 1191 |
|
| 1192 |
#ifdef IS_MPI |
| 1193 |
!! in MPI, we have to look up the unique IDs for each atom |
| 1194 |
unique_id_1 = AtomRowToGlobal(atom1) |
| 1195 |
unique_id_2 = AtomColToGlobal(atom2) |
| 1196 |
!! this situation should only arise in MPI simulations |
| 1197 |
if (unique_id_1 == unique_id_2) then |
| 1198 |
skip_it = .true. |
| 1199 |
return |
| 1200 |
end if |
| 1201 |
|
| 1202 |
!! this prevents us from doing the pair on multiple processors |
| 1203 |
if (unique_id_1 < unique_id_2) then |
| 1204 |
if (mod(unique_id_1 + unique_id_2,2) == 0) then |
| 1205 |
skip_it = .true. |
| 1206 |
return |
| 1207 |
endif |
| 1208 |
else |
| 1209 |
if (mod(unique_id_1 + unique_id_2,2) == 1) then |
| 1210 |
skip_it = .true. |
| 1211 |
return |
| 1212 |
endif |
| 1213 |
endif |
| 1214 |
#else |
| 1215 |
!! in the normal loop, the atom numbers are unique |
| 1216 |
unique_id_1 = atom1 |
| 1217 |
unique_id_2 = atom2 |
| 1218 |
#endif |
| 1219 |
|
| 1220 |
#ifdef IS_MPI |
| 1221 |
do i = 1, nSkipsForRowAtom(atom1) |
| 1222 |
if (skipsForRowAtom(atom1, i) .eq. unique_id_2) then |
| 1223 |
skip_it = .true. |
| 1224 |
return |
| 1225 |
endif |
| 1226 |
end do |
| 1227 |
#else |
| 1228 |
do i = 1, nSkipsForLocalAtom(atom1) |
| 1229 |
if (skipsForLocalAtom(atom1, i) .eq. unique_id_2) then |
| 1230 |
skip_it = .true. |
| 1231 |
return |
| 1232 |
endif |
| 1233 |
end do |
| 1234 |
#endif |
| 1235 |
|
| 1236 |
return |
| 1237 |
end function skipThisPair |
| 1238 |
|
| 1239 |
function getTopoDistance(atom1, atom2) result(topoDist) |
| 1240 |
integer, intent(in) :: atom1 |
| 1241 |
integer, intent(in) :: atom2 |
| 1242 |
integer :: topoDist |
| 1243 |
integer :: unique_id_2 |
| 1244 |
integer :: i |
| 1245 |
|
| 1246 |
#ifdef IS_MPI |
| 1247 |
unique_id_2 = AtomColToGlobal(atom2) |
| 1248 |
#else |
| 1249 |
unique_id_2 = atom2 |
| 1250 |
#endif |
| 1251 |
|
| 1252 |
! zero is default for unconnected (i.e. normal) pair interactions |
| 1253 |
|
| 1254 |
topoDist = 0 |
| 1255 |
|
| 1256 |
do i = 1, nTopoPairsForAtom(atom1) |
| 1257 |
if (toposForAtom(atom1, i) .eq. unique_id_2) then |
| 1258 |
topoDist = topoDistance(atom1, i) |
| 1259 |
return |
| 1260 |
endif |
| 1261 |
end do |
| 1262 |
|
| 1263 |
return |
| 1264 |
end function getTopoDistance |
| 1265 |
|
| 1266 |
function FF_UsesDirectionalAtoms() result(doesit) |
| 1267 |
logical :: doesit |
| 1268 |
doesit = FF_uses_DirectionalAtoms |
| 1269 |
end function FF_UsesDirectionalAtoms |
| 1270 |
|
| 1271 |
function FF_RequiresPrepairCalc() result(doesit) |
| 1272 |
logical :: doesit |
| 1273 |
doesit = FF_uses_EAM .or. FF_uses_SC |
| 1274 |
end function FF_RequiresPrepairCalc |
| 1275 |
|
| 1276 |
#ifdef PROFILE |
| 1277 |
function getforcetime() result(totalforcetime) |
| 1278 |
real(kind=dp) :: totalforcetime |
| 1279 |
totalforcetime = forcetime |
| 1280 |
end function getforcetime |
| 1281 |
#endif |
| 1282 |
|
| 1283 |
!! This cleans componets of force arrays belonging only to fortran |
| 1284 |
|
| 1285 |
subroutine add_stress_tensor(dpair, fpair, tau) |
| 1286 |
|
| 1287 |
real( kind = dp ), dimension(3), intent(in) :: dpair, fpair |
| 1288 |
real( kind = dp ), dimension(9), intent(inout) :: tau |
| 1289 |
|
| 1290 |
! because the d vector is the rj - ri vector, and |
| 1291 |
! because fx, fy, fz are the force on atom i, we need a |
| 1292 |
! negative sign here: |
| 1293 |
|
| 1294 |
tau(1) = tau(1) - dpair(1) * fpair(1) |
| 1295 |
tau(2) = tau(2) - dpair(1) * fpair(2) |
| 1296 |
tau(3) = tau(3) - dpair(1) * fpair(3) |
| 1297 |
tau(4) = tau(4) - dpair(2) * fpair(1) |
| 1298 |
tau(5) = tau(5) - dpair(2) * fpair(2) |
| 1299 |
tau(6) = tau(6) - dpair(2) * fpair(3) |
| 1300 |
tau(7) = tau(7) - dpair(3) * fpair(1) |
| 1301 |
tau(8) = tau(8) - dpair(3) * fpair(2) |
| 1302 |
tau(9) = tau(9) - dpair(3) * fpair(3) |
| 1303 |
|
| 1304 |
end subroutine add_stress_tensor |
| 1305 |
|
| 1306 |
end module doForces |