| 1 |
chuckv |
702 |
!! |
| 2 |
|
|
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
| 3 |
|
|
!! |
| 4 |
|
|
!! The University of Notre Dame grants you ("Licensee") a |
| 5 |
|
|
!! non-exclusive, royalty free, license to use, modify and |
| 6 |
|
|
!! redistribute this software in source and binary code form, provided |
| 7 |
|
|
!! that the following conditions are met: |
| 8 |
|
|
!! |
| 9 |
|
|
!! 1. Acknowledgement of the program authors must be made in any |
| 10 |
|
|
!! publication of scientific results based in part on use of the |
| 11 |
|
|
!! program. An acceptable form of acknowledgement is citation of |
| 12 |
|
|
!! the article in which the program was described (Matthew |
| 13 |
|
|
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
| 14 |
|
|
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
| 15 |
|
|
!! Parallel Simulation Engine for Molecular Dynamics," |
| 16 |
|
|
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
| 17 |
|
|
!! |
| 18 |
|
|
!! 2. Redistributions of source code must retain the above copyright |
| 19 |
|
|
!! notice, this list of conditions and the following disclaimer. |
| 20 |
|
|
!! |
| 21 |
|
|
!! 3. Redistributions in binary form must reproduce the above copyright |
| 22 |
|
|
!! notice, this list of conditions and the following disclaimer in the |
| 23 |
|
|
!! documentation and/or other materials provided with the |
| 24 |
|
|
!! distribution. |
| 25 |
|
|
!! |
| 26 |
|
|
!! This software is provided "AS IS," without a warranty of any |
| 27 |
|
|
!! kind. All express or implied conditions, representations and |
| 28 |
|
|
!! warranties, including any implied warranty of merchantability, |
| 29 |
|
|
!! fitness for a particular purpose or non-infringement, are hereby |
| 30 |
|
|
!! excluded. The University of Notre Dame and its licensors shall not |
| 31 |
|
|
!! be liable for any damages suffered by licensee as a result of |
| 32 |
|
|
!! using, modifying or distributing the software or its |
| 33 |
|
|
!! derivatives. In no event will the University of Notre Dame or its |
| 34 |
|
|
!! licensors be liable for any lost revenue, profit or data, or for |
| 35 |
|
|
!! direct, indirect, special, consequential, incidental or punitive |
| 36 |
|
|
!! damages, however caused and regardless of the theory of liability, |
| 37 |
|
|
!! arising out of the use of or inability to use software, even if the |
| 38 |
|
|
!! University of Notre Dame has been advised of the possibility of |
| 39 |
|
|
!! such damages. |
| 40 |
|
|
!! |
| 41 |
|
|
|
| 42 |
|
|
!! Impliments Sutton-Chen Metallic Potential |
| 43 |
|
|
!! See A.P.SUTTON and J.CHEN,PHIL MAG LETT 61,139-146,1990 |
| 44 |
|
|
|
| 45 |
|
|
|
| 46 |
|
|
module suttonchen |
| 47 |
|
|
use simulation |
| 48 |
|
|
use force_globals |
| 49 |
|
|
use status |
| 50 |
|
|
use atype_module |
| 51 |
|
|
use vector_class |
| 52 |
chuckv |
798 |
use fForceOptions |
| 53 |
chuckv |
702 |
#ifdef IS_MPI |
| 54 |
|
|
use mpiSimulation |
| 55 |
|
|
#endif |
| 56 |
|
|
implicit none |
| 57 |
|
|
PRIVATE |
| 58 |
|
|
#define __FORTRAN90 |
| 59 |
|
|
#include "UseTheForce/DarkSide/fInteractionMap.h" |
| 60 |
|
|
|
| 61 |
|
|
INTEGER, PARAMETER :: DP = selected_real_kind(15) |
| 62 |
|
|
|
| 63 |
|
|
logical, save :: SC_FF_initialized = .false. |
| 64 |
|
|
integer, save :: SC_Mixing_Policy |
| 65 |
|
|
real(kind = dp), save :: SC_rcut |
| 66 |
|
|
logical, save :: haveRcut = .false. |
| 67 |
chuckv |
728 |
logical, save :: haveMixingMap = .false. |
| 68 |
|
|
logical, save :: useGeometricDistanceMixing = .false. |
| 69 |
chuckv |
835 |
logical, save :: cleanArrays = .true. |
| 70 |
|
|
logical, save :: arraysAllocated = .false. |
| 71 |
chuckv |
702 |
|
| 72 |
chuckv |
728 |
|
| 73 |
chuckv |
702 |
character(len = statusMsgSize) :: errMesg |
| 74 |
chuckv |
714 |
integer :: sc_err |
| 75 |
chuckv |
702 |
|
| 76 |
|
|
character(len = 200) :: errMsg |
| 77 |
|
|
character(len=*), parameter :: RoutineName = "Sutton-Chen MODULE" |
| 78 |
|
|
!! Logical that determines if eam arrays should be zeroed |
| 79 |
|
|
logical :: nmflag = .false. |
| 80 |
|
|
|
| 81 |
|
|
|
| 82 |
|
|
type, private :: SCtype |
| 83 |
|
|
integer :: atid |
| 84 |
chuckv |
707 |
real(kind=dp) :: c |
| 85 |
chuckv |
702 |
real(kind=dp) :: m |
| 86 |
|
|
real(kind=dp) :: n |
| 87 |
|
|
real(kind=dp) :: alpha |
| 88 |
|
|
real(kind=dp) :: epsilon |
| 89 |
chuckv |
728 |
real(kind=dp) :: sc_rcut |
| 90 |
chuckv |
702 |
end type SCtype |
| 91 |
|
|
|
| 92 |
|
|
|
| 93 |
|
|
!! Arrays for derivatives used in force calculation |
| 94 |
chuckv |
713 |
real( kind = dp), dimension(:), allocatable :: rho |
| 95 |
chuckv |
702 |
real( kind = dp), dimension(:), allocatable :: frho |
| 96 |
|
|
real( kind = dp), dimension(:), allocatable :: dfrhodrho |
| 97 |
|
|
|
| 98 |
|
|
|
| 99 |
chuckv |
728 |
|
| 100 |
chuckv |
702 |
!! Arrays for MPI storage |
| 101 |
|
|
#ifdef IS_MPI |
| 102 |
|
|
real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_col |
| 103 |
|
|
real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_row |
| 104 |
|
|
real( kind = dp),save, dimension(:), allocatable :: frho_row |
| 105 |
|
|
real( kind = dp),save, dimension(:), allocatable :: frho_col |
| 106 |
|
|
real( kind = dp),save, dimension(:), allocatable :: rho_row |
| 107 |
|
|
real( kind = dp),save, dimension(:), allocatable :: rho_col |
| 108 |
|
|
real( kind = dp),save, dimension(:), allocatable :: rho_tmp |
| 109 |
|
|
#endif |
| 110 |
|
|
|
| 111 |
|
|
type, private :: SCTypeList |
| 112 |
|
|
integer :: nSCTypes = 0 |
| 113 |
|
|
integer :: currentSCtype = 0 |
| 114 |
|
|
|
| 115 |
|
|
type (SCtype), pointer :: SCtypes(:) => null() |
| 116 |
|
|
integer, pointer :: atidToSCtype(:) => null() |
| 117 |
|
|
end type SCTypeList |
| 118 |
|
|
|
| 119 |
|
|
type (SCTypeList), save :: SCList |
| 120 |
|
|
|
| 121 |
|
|
|
| 122 |
|
|
|
| 123 |
chuckv |
707 |
|
| 124 |
|
|
type :: MixParameters |
| 125 |
|
|
real(kind=DP) :: alpha |
| 126 |
|
|
real(kind=DP) :: epsilon |
| 127 |
chuckv |
728 |
real(kind=DP) :: m |
| 128 |
|
|
real(Kind=DP) :: n |
| 129 |
|
|
real(kind=DP) :: vpair_pot |
| 130 |
chuckv |
707 |
real(kind=dp) :: rCut |
| 131 |
|
|
logical :: rCutWasSet = .false. |
| 132 |
chuckv |
728 |
|
| 133 |
chuckv |
707 |
end type MixParameters |
| 134 |
|
|
|
| 135 |
|
|
type(MixParameters), dimension(:,:), allocatable :: MixingMap |
| 136 |
|
|
|
| 137 |
|
|
|
| 138 |
|
|
|
| 139 |
chuckv |
702 |
public :: setCutoffSC |
| 140 |
|
|
public :: do_SC_pair |
| 141 |
|
|
public :: newSCtype |
| 142 |
|
|
public :: calc_SC_prepair_rho |
| 143 |
chuckv |
735 |
public :: calc_SC_preforce_Frho |
| 144 |
chuckv |
702 |
public :: clean_SC |
| 145 |
|
|
public :: destroySCtypes |
| 146 |
|
|
public :: getSCCut |
| 147 |
chuckv |
728 |
! public :: setSCDefaultCutoff |
| 148 |
|
|
! public :: setSCUniformCutoff |
| 149 |
chuckv |
798 |
|
| 150 |
chuckv |
702 |
|
| 151 |
|
|
contains |
| 152 |
|
|
|
| 153 |
|
|
|
| 154 |
chuckv |
728 |
subroutine newSCtype(c_ident,c,m,n,alpha,epsilon,status) |
| 155 |
|
|
real (kind = dp ) :: c ! Density Scaling |
| 156 |
|
|
real (kind = dp ) :: m ! Density Exponent |
| 157 |
|
|
real (kind = dp ) :: n ! Pair Potential Exponent |
| 158 |
|
|
real (kind = dp ) :: alpha ! Length Scaling |
| 159 |
|
|
real (kind = dp ) :: epsilon ! Energy Scaling |
| 160 |
|
|
|
| 161 |
|
|
|
| 162 |
chuckv |
702 |
integer :: c_ident |
| 163 |
|
|
integer :: status |
| 164 |
|
|
|
| 165 |
chuckv |
713 |
integer :: nAtypes,nSCTypes,myATID |
| 166 |
chuckv |
702 |
integer :: maxVals |
| 167 |
|
|
integer :: alloc_stat |
| 168 |
|
|
integer :: current |
| 169 |
|
|
integer,pointer :: Matchlist(:) => null() |
| 170 |
|
|
|
| 171 |
|
|
status = 0 |
| 172 |
|
|
|
| 173 |
|
|
|
| 174 |
|
|
!! Assume that atypes has already been set and get the total number of types in atypes |
| 175 |
|
|
|
| 176 |
|
|
|
| 177 |
|
|
! check to see if this is the first time into |
| 178 |
chuckv |
728 |
if (.not.associated(SCList%SCTypes)) then |
| 179 |
chuckv |
831 |
call getMatchingElementList(atypes, "is_SC", .true., nSCtypes, MatchList) |
| 180 |
chuckv |
728 |
SCList%nSCtypes = nSCtypes |
| 181 |
|
|
allocate(SCList%SCTypes(nSCTypes)) |
| 182 |
chuckv |
702 |
nAtypes = getSize(atypes) |
| 183 |
chuckv |
728 |
allocate(SCList%atidToSCType(nAtypes)) |
| 184 |
chuckv |
702 |
end if |
| 185 |
|
|
|
| 186 |
chuckv |
728 |
SCList%currentSCType = SCList%currentSCType + 1 |
| 187 |
|
|
current = SCList%currentSCType |
| 188 |
chuckv |
702 |
|
| 189 |
|
|
myATID = getFirstMatchingElement(atypes, "c_ident", c_ident) |
| 190 |
chuckv |
728 |
SCList%atidToSCType(myATID) = current |
| 191 |
chuckv |
702 |
|
| 192 |
|
|
|
| 193 |
|
|
|
| 194 |
chuckv |
728 |
SCList%SCTypes(current)%atid = c_ident |
| 195 |
|
|
SCList%SCTypes(current)%alpha = alpha |
| 196 |
|
|
SCList%SCTypes(current)%c = c |
| 197 |
|
|
SCList%SCTypes(current)%m = m |
| 198 |
|
|
SCList%SCTypes(current)%n = n |
| 199 |
|
|
SCList%SCTypes(current)%epsilon = epsilon |
| 200 |
chuckv |
713 |
end subroutine newSCtype |
| 201 |
chuckv |
702 |
|
| 202 |
chuckv |
713 |
|
| 203 |
|
|
subroutine destroySCTypes() |
| 204 |
chuckv |
728 |
if (associated(SCList%SCtypes)) then |
| 205 |
|
|
deallocate(SCList%SCTypes) |
| 206 |
|
|
SCList%SCTypes=>null() |
| 207 |
|
|
end if |
| 208 |
|
|
if (associated(SCList%atidToSCtype)) then |
| 209 |
|
|
deallocate(SCList%atidToSCtype) |
| 210 |
|
|
SCList%atidToSCtype=>null() |
| 211 |
|
|
end if |
| 212 |
chuckv |
702 |
|
| 213 |
|
|
|
| 214 |
chuckv |
713 |
end subroutine destroySCTypes |
| 215 |
chuckv |
702 |
|
| 216 |
|
|
|
| 217 |
chuckv |
713 |
|
| 218 |
|
|
function getSCCut(atomID) result(cutValue) |
| 219 |
chuckv |
702 |
integer, intent(in) :: atomID |
| 220 |
chuckv |
728 |
integer :: scID |
| 221 |
chuckv |
702 |
real(kind=dp) :: cutValue |
| 222 |
|
|
|
| 223 |
chuckv |
728 |
scID = SCList%atidToSCType(atomID) |
| 224 |
chuckv |
831 |
cutValue = 2.0_dp * SCList%SCTypes(scID)%alpha |
| 225 |
chuckv |
713 |
end function getSCCut |
| 226 |
chuckv |
702 |
|
| 227 |
|
|
|
| 228 |
|
|
|
| 229 |
chuckv |
713 |
|
| 230 |
|
|
subroutine createMixingMap() |
| 231 |
|
|
integer :: nSCtypes, i, j |
| 232 |
|
|
real ( kind = dp ) :: e1, e2,m1,m2,alpha1,alpha2,n1,n2 |
| 233 |
|
|
real ( kind = dp ) :: rcut6, tp6, tp12 |
| 234 |
|
|
logical :: isSoftCore1, isSoftCore2, doShift |
| 235 |
|
|
|
| 236 |
|
|
if (SCList%currentSCtype == 0) then |
| 237 |
|
|
call handleError("SuttonChen", "No members in SCMap") |
| 238 |
chuckv |
702 |
return |
| 239 |
|
|
end if |
| 240 |
|
|
|
| 241 |
chuckv |
713 |
nSCtypes = SCList%nSCtypes |
| 242 |
chuckv |
702 |
|
| 243 |
chuckv |
713 |
if (.not. allocated(MixingMap)) then |
| 244 |
|
|
allocate(MixingMap(nSCtypes, nSCtypes)) |
| 245 |
|
|
endif |
| 246 |
chuckv |
798 |
useGeometricDistanceMixing = usesGeometricDistanceMixing() |
| 247 |
chuckv |
713 |
do i = 1, nSCtypes |
| 248 |
chuckv |
702 |
|
| 249 |
chuckv |
713 |
e1 = SCList%SCtypes(i)%epsilon |
| 250 |
|
|
m1 = SCList%SCtypes(i)%m |
| 251 |
|
|
n1 = SCList%SCtypes(i)%n |
| 252 |
|
|
alpha1 = SCList%SCtypes(i)%alpha |
| 253 |
chuckv |
702 |
|
| 254 |
chuckv |
713 |
do j = i, nSCtypes |
| 255 |
|
|
|
| 256 |
chuckv |
702 |
|
| 257 |
chuckv |
713 |
e2 = SCList%SCtypes(j)%epsilon |
| 258 |
|
|
m2 = SCList%SCtypes(j)%m |
| 259 |
|
|
n2 = SCList%SCtypes(j)%n |
| 260 |
|
|
alpha2 = SCList%SCtypes(j)%alpha |
| 261 |
chuckv |
702 |
|
| 262 |
chuckv |
728 |
if (useGeometricDistanceMixing) then |
| 263 |
|
|
MixingMap(i,j)%alpha = sqrt(alpha1 * alpha2) !SC formulation |
| 264 |
|
|
else |
| 265 |
|
|
MixingMap(i,j)%alpha = 0.5_dp * (alpha1 + alpha2) ! Goddard formulation |
| 266 |
|
|
endif |
| 267 |
|
|
|
| 268 |
|
|
MixingMap(i,j)%epsilon = sqrt(e1 * e2) |
| 269 |
chuckv |
713 |
MixingMap(i,j)%m = 0.5_dp*(m1+m2) |
| 270 |
|
|
MixingMap(i,j)%n = 0.5_dp*(n1+n2) |
| 271 |
|
|
MixingMap(i,j)%alpha = 0.5_dp*(alpha1+alpha2) |
| 272 |
|
|
MixingMap(i,j)%rcut = 2.0_dp *MixingMap(i,j)%alpha |
| 273 |
chuckv |
728 |
MixingMap(i,j)%vpair_pot = MixingMap(i,j)%epsilon* & |
| 274 |
|
|
(MixingMap(i,j)%alpha/MixingMap(i,j)%rcut)**MixingMap(i,j)%n |
| 275 |
|
|
if (i.ne.j) then |
| 276 |
|
|
MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon |
| 277 |
|
|
MixingMap(j,i)%m = MixingMap(i,j)%m |
| 278 |
|
|
MixingMap(j,i)%n = MixingMap(i,j)%n |
| 279 |
|
|
MixingMap(j,i)%alpha = MixingMap(i,j)%alpha |
| 280 |
|
|
MixingMap(j,i)%rcut = MixingMap(i,j)%rcut |
| 281 |
|
|
MixingMap(j,i)%vpair_pot = MixingMap(i,j)%vpair_pot |
| 282 |
|
|
endif |
| 283 |
chuckv |
713 |
enddo |
| 284 |
chuckv |
702 |
enddo |
| 285 |
chuckv |
713 |
|
| 286 |
|
|
haveMixingMap = .true. |
| 287 |
|
|
|
| 288 |
|
|
end subroutine createMixingMap |
| 289 |
|
|
|
| 290 |
chuckv |
702 |
|
| 291 |
|
|
!! routine checks to see if array is allocated, deallocates array if allocated |
| 292 |
|
|
!! and then creates the array to the required size |
| 293 |
chuckv |
835 |
subroutine allocateSC() |
| 294 |
|
|
integer :: status |
| 295 |
chuckv |
702 |
|
| 296 |
|
|
#ifdef IS_MPI |
| 297 |
|
|
integer :: nAtomsInRow |
| 298 |
|
|
integer :: nAtomsInCol |
| 299 |
|
|
#endif |
| 300 |
|
|
integer :: alloc_stat |
| 301 |
|
|
|
| 302 |
chuckv |
835 |
|
| 303 |
chuckv |
702 |
status = 0 |
| 304 |
|
|
#ifdef IS_MPI |
| 305 |
|
|
nAtomsInRow = getNatomsInRow(plan_atom_row) |
| 306 |
|
|
nAtomsInCol = getNatomsInCol(plan_atom_col) |
| 307 |
|
|
#endif |
| 308 |
|
|
|
| 309 |
chuckv |
713 |
|
| 310 |
|
|
|
| 311 |
chuckv |
702 |
if (allocated(frho)) deallocate(frho) |
| 312 |
|
|
allocate(frho(nlocal),stat=alloc_stat) |
| 313 |
chuckv |
713 |
if (alloc_stat /= 0) then |
| 314 |
chuckv |
702 |
status = -1 |
| 315 |
|
|
end if |
| 316 |
chuckv |
713 |
|
| 317 |
chuckv |
702 |
if (allocated(rho)) deallocate(rho) |
| 318 |
|
|
allocate(rho(nlocal),stat=alloc_stat) |
| 319 |
|
|
if (alloc_stat /= 0) then |
| 320 |
|
|
status = -1 |
| 321 |
|
|
end if |
| 322 |
|
|
|
| 323 |
|
|
if (allocated(dfrhodrho)) deallocate(dfrhodrho) |
| 324 |
|
|
allocate(dfrhodrho(nlocal),stat=alloc_stat) |
| 325 |
|
|
if (alloc_stat /= 0) then |
| 326 |
|
|
status = -1 |
| 327 |
|
|
end if |
| 328 |
|
|
|
| 329 |
|
|
#ifdef IS_MPI |
| 330 |
|
|
|
| 331 |
|
|
if (allocated(rho_tmp)) deallocate(rho_tmp) |
| 332 |
|
|
allocate(rho_tmp(nlocal),stat=alloc_stat) |
| 333 |
|
|
if (alloc_stat /= 0) then |
| 334 |
|
|
status = -1 |
| 335 |
|
|
end if |
| 336 |
|
|
|
| 337 |
|
|
|
| 338 |
|
|
if (allocated(frho_row)) deallocate(frho_row) |
| 339 |
|
|
allocate(frho_row(nAtomsInRow),stat=alloc_stat) |
| 340 |
|
|
if (alloc_stat /= 0) then |
| 341 |
|
|
status = -1 |
| 342 |
|
|
end if |
| 343 |
|
|
if (allocated(rho_row)) deallocate(rho_row) |
| 344 |
|
|
allocate(rho_row(nAtomsInRow),stat=alloc_stat) |
| 345 |
|
|
if (alloc_stat /= 0) then |
| 346 |
|
|
status = -1 |
| 347 |
|
|
end if |
| 348 |
|
|
if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row) |
| 349 |
|
|
allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat) |
| 350 |
|
|
if (alloc_stat /= 0) then |
| 351 |
|
|
status = -1 |
| 352 |
|
|
end if |
| 353 |
|
|
|
| 354 |
|
|
|
| 355 |
|
|
! Now do column arrays |
| 356 |
|
|
|
| 357 |
|
|
if (allocated(frho_col)) deallocate(frho_col) |
| 358 |
|
|
allocate(frho_col(nAtomsInCol),stat=alloc_stat) |
| 359 |
|
|
if (alloc_stat /= 0) then |
| 360 |
|
|
status = -1 |
| 361 |
|
|
end if |
| 362 |
|
|
if (allocated(rho_col)) deallocate(rho_col) |
| 363 |
|
|
allocate(rho_col(nAtomsInCol),stat=alloc_stat) |
| 364 |
|
|
if (alloc_stat /= 0) then |
| 365 |
|
|
status = -1 |
| 366 |
|
|
end if |
| 367 |
|
|
if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col) |
| 368 |
|
|
allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat) |
| 369 |
|
|
if (alloc_stat /= 0) then |
| 370 |
|
|
status = -1 |
| 371 |
|
|
end if |
| 372 |
|
|
|
| 373 |
|
|
#endif |
| 374 |
chuckv |
835 |
if (status == -1) then |
| 375 |
|
|
call handleError("SuttonChen:allocateSC","Error in allocating SC arrays") |
| 376 |
|
|
end if |
| 377 |
|
|
arraysAllocated = .true. |
| 378 |
chuckv |
713 |
end subroutine allocateSC |
| 379 |
chuckv |
702 |
|
| 380 |
|
|
!! C sets rcut to be the largest cutoff of any atype |
| 381 |
|
|
!! present in this simulation. Doesn't include all atypes |
| 382 |
|
|
!! sim knows about, just those in the simulation. |
| 383 |
chuckv |
713 |
subroutine setCutoffSC(rcut, status) |
| 384 |
chuckv |
702 |
real(kind=dp) :: rcut |
| 385 |
|
|
integer :: status |
| 386 |
|
|
status = 0 |
| 387 |
|
|
|
| 388 |
chuckv |
713 |
SC_rcut = rcut |
| 389 |
chuckv |
702 |
|
| 390 |
chuckv |
713 |
end subroutine setCutoffSC |
| 391 |
chuckv |
702 |
|
| 392 |
chuckv |
835 |
!! This array allocates module arrays if needed and builds mixing map. |
| 393 |
chuckv |
728 |
subroutine clean_SC() |
| 394 |
chuckv |
835 |
if (.not.arraysAllocated) call allocateSC() |
| 395 |
|
|
if (.not.haveMixingMap) call createMixingMap() |
| 396 |
chuckv |
702 |
! clean non-IS_MPI first |
| 397 |
|
|
frho = 0.0_dp |
| 398 |
|
|
rho = 0.0_dp |
| 399 |
|
|
dfrhodrho = 0.0_dp |
| 400 |
|
|
! clean MPI if needed |
| 401 |
|
|
#ifdef IS_MPI |
| 402 |
|
|
frho_row = 0.0_dp |
| 403 |
|
|
frho_col = 0.0_dp |
| 404 |
|
|
rho_row = 0.0_dp |
| 405 |
|
|
rho_col = 0.0_dp |
| 406 |
|
|
rho_tmp = 0.0_dp |
| 407 |
|
|
dfrhodrho_row = 0.0_dp |
| 408 |
|
|
dfrhodrho_col = 0.0_dp |
| 409 |
|
|
#endif |
| 410 |
chuckv |
728 |
end subroutine clean_SC |
| 411 |
chuckv |
702 |
|
| 412 |
|
|
|
| 413 |
|
|
|
| 414 |
|
|
!! Calculates rho_r |
| 415 |
gezelter |
762 |
subroutine calc_sc_prepair_rho(atom1,atom2,d,r,rijsq, rcut) |
| 416 |
chuckv |
702 |
integer :: atom1,atom2 |
| 417 |
|
|
real(kind = dp), dimension(3) :: d |
| 418 |
|
|
real(kind = dp), intent(inout) :: r |
| 419 |
gezelter |
762 |
real(kind = dp), intent(inout) :: rijsq, rcut |
| 420 |
chuckv |
702 |
! value of electron density rho do to atom i at atom j |
| 421 |
|
|
real(kind = dp) :: rho_i_at_j |
| 422 |
|
|
! value of electron density rho do to atom j at atom i |
| 423 |
|
|
real(kind = dp) :: rho_j_at_i |
| 424 |
|
|
|
| 425 |
|
|
! we don't use the derivatives, dummy variables |
| 426 |
chuckv |
728 |
real( kind = dp) :: drho |
| 427 |
chuckv |
714 |
integer :: sc_err |
| 428 |
chuckv |
702 |
|
| 429 |
|
|
integer :: atid1,atid2 ! Global atid |
| 430 |
|
|
integer :: myid_atom1 ! EAM atid |
| 431 |
|
|
integer :: myid_atom2 |
| 432 |
|
|
|
| 433 |
|
|
|
| 434 |
|
|
! check to see if we need to be cleaned at the start of a force loop |
| 435 |
|
|
|
| 436 |
chuckv |
835 |
if (cleanArrays) call clean_SC() |
| 437 |
|
|
cleanArrays = .false. |
| 438 |
chuckv |
702 |
|
| 439 |
|
|
#ifdef IS_MPI |
| 440 |
|
|
Atid1 = Atid_row(Atom1) |
| 441 |
|
|
Atid2 = Atid_col(Atom2) |
| 442 |
|
|
#else |
| 443 |
|
|
Atid1 = Atid(Atom1) |
| 444 |
|
|
Atid2 = Atid(Atom2) |
| 445 |
|
|
#endif |
| 446 |
|
|
|
| 447 |
chuckv |
728 |
Myid_atom1 = SCList%atidtoSCtype(Atid1) |
| 448 |
|
|
Myid_atom2 = SCList%atidtoSCtype(Atid2) |
| 449 |
chuckv |
702 |
|
| 450 |
|
|
|
| 451 |
|
|
|
| 452 |
chuckv |
713 |
rho_i_at_j = (MixingMap(Myid_atom1,Myid_atom2)%alpha/r)& |
| 453 |
|
|
**MixingMap(Myid_atom1,Myid_atom2)%m |
| 454 |
|
|
rho_j_at_i = rho_i_at_j |
| 455 |
chuckv |
702 |
|
| 456 |
|
|
#ifdef IS_MPI |
| 457 |
|
|
rho_col(atom2) = rho_col(atom2) + rho_i_at_j |
| 458 |
chuckv |
713 |
rho_row(atom1) = rho_row(atom1) + rho_j_at_i |
| 459 |
chuckv |
702 |
#else |
| 460 |
|
|
rho(atom2) = rho(atom2) + rho_i_at_j |
| 461 |
|
|
rho(atom1) = rho(atom1) + rho_j_at_i |
| 462 |
|
|
#endif |
| 463 |
|
|
|
| 464 |
chuckv |
713 |
end subroutine calc_sc_prepair_rho |
| 465 |
chuckv |
702 |
|
| 466 |
|
|
|
| 467 |
chuckv |
714 |
!! Calculate the rho_a for all local atoms |
| 468 |
chuckv |
728 |
subroutine calc_sc_preforce_Frho(nlocal,pot) |
| 469 |
chuckv |
702 |
integer :: nlocal |
| 470 |
|
|
real(kind=dp) :: pot |
| 471 |
|
|
integer :: i,j |
| 472 |
|
|
integer :: atom |
| 473 |
|
|
real(kind=dp) :: U,U1,U2 |
| 474 |
|
|
integer :: atype1 |
| 475 |
chuckv |
728 |
integer :: atid1 |
| 476 |
|
|
integer :: myid |
| 477 |
chuckv |
702 |
|
| 478 |
|
|
|
| 479 |
|
|
!! Scatter the electron density from pre-pair calculation back to local atoms |
| 480 |
|
|
#ifdef IS_MPI |
| 481 |
chuckv |
714 |
call scatter(rho_row,rho,plan_atom_row,sc_err) |
| 482 |
|
|
if (sc_err /= 0 ) then |
| 483 |
chuckv |
702 |
write(errMsg,*) " Error scattering rho_row into rho" |
| 484 |
|
|
call handleError(RoutineName,errMesg) |
| 485 |
|
|
endif |
| 486 |
chuckv |
714 |
call scatter(rho_col,rho_tmp,plan_atom_col,sc_err) |
| 487 |
|
|
if (sc_err /= 0 ) then |
| 488 |
chuckv |
702 |
write(errMsg,*) " Error scattering rho_col into rho" |
| 489 |
|
|
call handleError(RoutineName,errMesg) |
| 490 |
chuckv |
713 |
|
| 491 |
chuckv |
702 |
endif |
| 492 |
|
|
|
| 493 |
|
|
rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal) |
| 494 |
|
|
#endif |
| 495 |
|
|
|
| 496 |
|
|
|
| 497 |
|
|
|
| 498 |
chuckv |
713 |
!! Calculate F(rho) and derivative |
| 499 |
chuckv |
702 |
do atom = 1, nlocal |
| 500 |
chuckv |
728 |
Myid = SCList%atidtoSctype(Atid(atom)) |
| 501 |
|
|
frho(atom) = -SCList%SCTypes(Myid)%c * & |
| 502 |
|
|
SCList%SCTypes(Myid)%epsilon * sqrt(rho(i)) |
| 503 |
|
|
|
| 504 |
chuckv |
713 |
dfrhodrho(atom) = 0.5_dp*frho(atom)/rho(atom) |
| 505 |
chuckv |
702 |
pot = pot + u |
| 506 |
|
|
enddo |
| 507 |
|
|
|
| 508 |
chuckv |
731 |
#ifdef IS_MPI |
| 509 |
chuckv |
702 |
!! communicate f(rho) and derivatives back into row and column arrays |
| 510 |
chuckv |
714 |
call gather(frho,frho_row,plan_atom_row, sc_err) |
| 511 |
|
|
if (sc_err /= 0) then |
| 512 |
chuckv |
702 |
call handleError("cal_eam_forces()","MPI gather frho_row failure") |
| 513 |
|
|
endif |
| 514 |
chuckv |
714 |
call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err) |
| 515 |
|
|
if (sc_err /= 0) then |
| 516 |
chuckv |
702 |
call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure") |
| 517 |
|
|
endif |
| 518 |
chuckv |
714 |
call gather(frho,frho_col,plan_atom_col, sc_err) |
| 519 |
|
|
if (sc_err /= 0) then |
| 520 |
chuckv |
702 |
call handleError("cal_eam_forces()","MPI gather frho_col failure") |
| 521 |
|
|
endif |
| 522 |
chuckv |
714 |
call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err) |
| 523 |
|
|
if (sc_err /= 0) then |
| 524 |
chuckv |
702 |
call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure") |
| 525 |
|
|
endif |
| 526 |
|
|
#endif |
| 527 |
chuckv |
713 |
|
| 528 |
|
|
|
| 529 |
chuckv |
728 |
end subroutine calc_sc_preforce_Frho |
| 530 |
chuckv |
702 |
|
| 531 |
|
|
|
| 532 |
chuckv |
713 |
!! Does Sutton-Chen pairwise Force calculation. |
| 533 |
gezelter |
762 |
subroutine do_sc_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, & |
| 534 |
chuckv |
702 |
pot, f, do_pot) |
| 535 |
|
|
!Arguments |
| 536 |
|
|
integer, intent(in) :: atom1, atom2 |
| 537 |
gezelter |
762 |
real( kind = dp ), intent(in) :: rij, r2, rcut |
| 538 |
chuckv |
702 |
real( kind = dp ) :: pot, sw, vpair |
| 539 |
|
|
real( kind = dp ), dimension(3,nLocal) :: f |
| 540 |
|
|
real( kind = dp ), intent(in), dimension(3) :: d |
| 541 |
|
|
real( kind = dp ), intent(inout), dimension(3) :: fpair |
| 542 |
|
|
|
| 543 |
|
|
logical, intent(in) :: do_pot |
| 544 |
|
|
|
| 545 |
|
|
real( kind = dp ) :: drdx,drdy,drdz |
| 546 |
chuckv |
728 |
real( kind = dp ) :: dvpdr |
| 547 |
|
|
real( kind = dp ) :: drhodr |
| 548 |
chuckv |
702 |
real( kind = dp ) :: dudr |
| 549 |
chuckv |
714 |
real( kind = dp ) :: rcij |
| 550 |
chuckv |
702 |
real( kind = dp ) :: drhoidr,drhojdr |
| 551 |
|
|
real( kind = dp ) :: Fx,Fy,Fz |
| 552 |
|
|
real( kind = dp ) :: r,d2pha,phb,d2phb |
| 553 |
chuckv |
728 |
real( kind = dp ) :: pot_temp,vptmp |
| 554 |
|
|
real( kind = dp ) :: epsilonij,aij,nij,mij,vcij |
| 555 |
chuckv |
702 |
integer :: id1,id2 |
| 556 |
|
|
integer :: mytype_atom1 |
| 557 |
|
|
integer :: mytype_atom2 |
| 558 |
|
|
integer :: atid1,atid2 |
| 559 |
|
|
!Local Variables |
| 560 |
|
|
|
| 561 |
|
|
! write(*,*) "Frho: ", Frho(atom1) |
| 562 |
chuckv |
835 |
|
| 563 |
|
|
cleanArrays = .true. |
| 564 |
chuckv |
702 |
|
| 565 |
|
|
dvpdr = 0.0E0_DP |
| 566 |
|
|
|
| 567 |
|
|
|
| 568 |
|
|
#ifdef IS_MPI |
| 569 |
|
|
atid1 = atid_row(atom1) |
| 570 |
|
|
atid2 = atid_col(atom2) |
| 571 |
|
|
#else |
| 572 |
|
|
atid1 = atid(atom1) |
| 573 |
|
|
atid2 = atid(atom2) |
| 574 |
|
|
#endif |
| 575 |
|
|
|
| 576 |
chuckv |
728 |
mytype_atom1 = SCList%atidToSCType(atid1) |
| 577 |
|
|
mytype_atom2 = SCList%atidTOSCType(atid2) |
| 578 |
chuckv |
702 |
|
| 579 |
|
|
drdx = d(1)/rij |
| 580 |
|
|
drdy = d(2)/rij |
| 581 |
|
|
drdz = d(3)/rij |
| 582 |
|
|
|
| 583 |
chuckv |
714 |
|
| 584 |
|
|
epsilonij = MixingMap(mytype_atom1,mytype_atom2)%epsilon |
| 585 |
|
|
aij = MixingMap(mytype_atom1,mytype_atom2)%alpha |
| 586 |
|
|
nij = MixingMap(mytype_atom1,mytype_atom2)%n |
| 587 |
|
|
mij = MixingMap(mytype_atom1,mytype_atom2)%m |
| 588 |
|
|
vcij = MixingMap(mytype_atom1,mytype_atom2)%vpair_pot |
| 589 |
chuckv |
702 |
|
| 590 |
chuckv |
728 |
vptmp = epsilonij*((aij/r)**nij) |
| 591 |
chuckv |
702 |
|
| 592 |
|
|
|
| 593 |
chuckv |
714 |
dvpdr = -nij*vptmp/r |
| 594 |
|
|
drhodr = -mij*((aij/r)**mij)/r |
| 595 |
chuckv |
728 |
|
| 596 |
chuckv |
714 |
|
| 597 |
chuckv |
728 |
dudr = drhodr*(dfrhodrho(atom1)+dfrhodrho(atom2)) & |
| 598 |
chuckv |
714 |
+ dvpdr |
| 599 |
chuckv |
728 |
|
| 600 |
|
|
pot_temp = vptmp + vcij |
| 601 |
chuckv |
714 |
|
| 602 |
chuckv |
702 |
|
| 603 |
|
|
#ifdef IS_MPI |
| 604 |
chuckv |
714 |
dudr = drhodr*(dfrhodrho_row(atom1)+dfrhodrho_col(atom2)) & |
| 605 |
chuckv |
702 |
+ dvpdr |
| 606 |
|
|
|
| 607 |
|
|
#else |
| 608 |
chuckv |
714 |
dudr = drhodr*(dfrhodrho(atom1)+dfrhodrho(atom2)) & |
| 609 |
chuckv |
702 |
+ dvpdr |
| 610 |
|
|
#endif |
| 611 |
|
|
|
| 612 |
chuckv |
714 |
|
| 613 |
chuckv |
702 |
fx = dudr * drdx |
| 614 |
|
|
fy = dudr * drdy |
| 615 |
|
|
fz = dudr * drdz |
| 616 |
|
|
|
| 617 |
|
|
|
| 618 |
|
|
#ifdef IS_MPI |
| 619 |
|
|
if (do_pot) then |
| 620 |
chuckv |
728 |
pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + (pot_temp)*0.5 |
| 621 |
|
|
pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + (pot_temp)*0.5 |
| 622 |
chuckv |
702 |
end if |
| 623 |
|
|
|
| 624 |
|
|
f_Row(1,atom1) = f_Row(1,atom1) + fx |
| 625 |
|
|
f_Row(2,atom1) = f_Row(2,atom1) + fy |
| 626 |
|
|
f_Row(3,atom1) = f_Row(3,atom1) + fz |
| 627 |
|
|
|
| 628 |
|
|
f_Col(1,atom2) = f_Col(1,atom2) - fx |
| 629 |
|
|
f_Col(2,atom2) = f_Col(2,atom2) - fy |
| 630 |
|
|
f_Col(3,atom2) = f_Col(3,atom2) - fz |
| 631 |
|
|
#else |
| 632 |
|
|
|
| 633 |
|
|
if(do_pot) then |
| 634 |
chuckv |
728 |
pot = pot + pot_temp |
| 635 |
chuckv |
702 |
end if |
| 636 |
|
|
|
| 637 |
|
|
f(1,atom1) = f(1,atom1) + fx |
| 638 |
|
|
f(2,atom1) = f(2,atom1) + fy |
| 639 |
|
|
f(3,atom1) = f(3,atom1) + fz |
| 640 |
|
|
|
| 641 |
|
|
f(1,atom2) = f(1,atom2) - fx |
| 642 |
|
|
f(2,atom2) = f(2,atom2) - fy |
| 643 |
|
|
f(3,atom2) = f(3,atom2) - fz |
| 644 |
|
|
#endif |
| 645 |
|
|
|
| 646 |
chuckv |
728 |
|
| 647 |
chuckv |
702 |
#ifdef IS_MPI |
| 648 |
|
|
id1 = AtomRowToGlobal(atom1) |
| 649 |
|
|
id2 = AtomColToGlobal(atom2) |
| 650 |
|
|
#else |
| 651 |
|
|
id1 = atom1 |
| 652 |
|
|
id2 = atom2 |
| 653 |
|
|
#endif |
| 654 |
|
|
|
| 655 |
|
|
if (molMembershipList(id1) .ne. molMembershipList(id2)) then |
| 656 |
|
|
|
| 657 |
|
|
fpair(1) = fpair(1) + fx |
| 658 |
|
|
fpair(2) = fpair(2) + fy |
| 659 |
|
|
fpair(3) = fpair(3) + fz |
| 660 |
|
|
|
| 661 |
|
|
endif |
| 662 |
|
|
|
| 663 |
|
|
|
| 664 |
chuckv |
728 |
end subroutine do_sc_pair |
| 665 |
chuckv |
702 |
|
| 666 |
|
|
|
| 667 |
|
|
|
| 668 |
chuckv |
713 |
end module suttonchen |