| 1 | gezelter | 246 | !! | 
| 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 | gezelter | 115 | module eam | 
| 43 |  |  | use simulation | 
| 44 |  |  | use force_globals | 
| 45 |  |  | use status | 
| 46 |  |  | use atype_module | 
| 47 | chrisfen | 578 | use vector_class | 
| 48 | gezelter | 938 | use interpolation | 
| 49 | gezelter | 115 | #ifdef IS_MPI | 
| 50 |  |  | use mpiSimulation | 
| 51 |  |  | #endif | 
| 52 |  |  | implicit none | 
| 53 |  |  | PRIVATE | 
| 54 | chuckv | 656 | #define __FORTRAN90 | 
| 55 |  |  | #include "UseTheForce/DarkSide/fInteractionMap.h" | 
| 56 | gezelter | 115 |  | 
| 57 | chuckv | 131 | INTEGER, PARAMETER :: DP = selected_real_kind(15) | 
| 58 | gezelter | 115 |  | 
| 59 |  |  | logical, save :: EAM_FF_initialized = .false. | 
| 60 |  |  | integer, save :: EAM_Mixing_Policy | 
| 61 |  |  | real(kind = dp), save :: EAM_rcut | 
| 62 |  |  | logical, save :: haveRcut = .false. | 
| 63 |  |  |  | 
| 64 |  |  | character(len = statusMsgSize) :: errMesg | 
| 65 |  |  | integer :: eam_err | 
| 66 |  |  |  | 
| 67 |  |  | character(len = 200) :: errMsg | 
| 68 |  |  | character(len=*), parameter :: RoutineName =  "EAM MODULE" | 
| 69 | gezelter | 507 | !! Logical that determines if eam arrays should be zeroed | 
| 70 | gezelter | 115 | logical :: cleanme = .true. | 
| 71 |  |  |  | 
| 72 |  |  | type, private :: EAMtype | 
| 73 |  |  | integer           :: eam_atype | 
| 74 |  |  | real( kind = DP ) :: eam_lattice | 
| 75 |  |  | real( kind = DP ) :: eam_rcut | 
| 76 |  |  | integer           :: eam_atype_map | 
| 77 | gezelter | 938 | type(cubicSpline) :: rho | 
| 78 |  |  | type(cubicSpline) :: Z | 
| 79 |  |  | type(cubicSpline) :: F | 
| 80 |  |  | type(cubicSpline) :: phi | 
| 81 | gezelter | 115 | end type EAMtype | 
| 82 |  |  |  | 
| 83 |  |  | !! Arrays for derivatives used in force calculation | 
| 84 |  |  | real( kind = dp), dimension(:), allocatable :: frho | 
| 85 |  |  | real( kind = dp), dimension(:), allocatable :: rho | 
| 86 |  |  | real( kind = dp), dimension(:), allocatable :: dfrhodrho | 
| 87 |  |  |  | 
| 88 | gezelter | 507 | !! Arrays for MPI storage | 
| 89 | gezelter | 115 | #ifdef IS_MPI | 
| 90 |  |  | real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_col | 
| 91 |  |  | real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_row | 
| 92 |  |  | real( kind = dp),save, dimension(:), allocatable :: frho_row | 
| 93 |  |  | real( kind = dp),save, dimension(:), allocatable :: frho_col | 
| 94 |  |  | real( kind = dp),save, dimension(:), allocatable :: rho_row | 
| 95 |  |  | real( kind = dp),save, dimension(:), allocatable :: rho_col | 
| 96 |  |  | real( kind = dp),save, dimension(:), allocatable :: rho_tmp | 
| 97 |  |  | #endif | 
| 98 |  |  |  | 
| 99 |  |  | type, private :: EAMTypeList | 
| 100 |  |  | integer           :: n_eam_types = 0 | 
| 101 |  |  | integer           :: currentAddition = 0 | 
| 102 | gezelter | 507 |  | 
| 103 | gezelter | 115 | type (EAMtype), pointer  :: EAMParams(:) => null() | 
| 104 | chuckv | 577 | integer, pointer         :: atidToEAMType(:) => null() | 
| 105 | gezelter | 115 | end type EAMTypeList | 
| 106 |  |  |  | 
| 107 |  |  | type (eamTypeList), save :: EAMList | 
| 108 |  |  |  | 
| 109 |  |  | !! standard eam stuff | 
| 110 |  |  |  | 
| 111 |  |  |  | 
| 112 |  |  | public :: init_EAM_FF | 
| 113 |  |  | public :: setCutoffEAM | 
| 114 |  |  | public :: do_eam_pair | 
| 115 |  |  | public :: newEAMtype | 
| 116 |  |  | public :: calc_eam_prepair_rho | 
| 117 |  |  | public :: calc_eam_preforce_Frho | 
| 118 |  |  | public :: clean_EAM | 
| 119 | chuckv | 472 | public :: destroyEAMTypes | 
| 120 | chrisfen | 578 | public :: getEAMCut | 
| 121 | chrisfen | 943 | public :: lookupEAMSpline | 
| 122 |  |  | public :: lookupEAMSpline1d | 
| 123 | gezelter | 115 |  | 
| 124 |  |  | contains | 
| 125 |  |  |  | 
| 126 |  |  | subroutine newEAMtype(lattice_constant,eam_nrho,eam_drho,eam_nr,& | 
| 127 | gezelter | 938 | eam_dr,rcut,eam_Z_r,eam_rho_r,eam_F_rho, c_ident, status) | 
| 128 | gezelter | 115 | real (kind = dp )                      :: lattice_constant | 
| 129 |  |  | integer                                :: eam_nrho | 
| 130 |  |  | real (kind = dp )                      :: eam_drho | 
| 131 |  |  | integer                                :: eam_nr | 
| 132 |  |  | real (kind = dp )                      :: eam_dr | 
| 133 |  |  | real (kind = dp )                      :: rcut | 
| 134 | gezelter | 938 | real (kind = dp ), dimension(eam_nr)   :: eam_Z_r, rvals | 
| 135 |  |  | real (kind = dp ), dimension(eam_nr)   :: eam_rho_r, eam_phi_r | 
| 136 |  |  | real (kind = dp ), dimension(eam_nrho) :: eam_F_rho, rhovals | 
| 137 | chrisfen | 579 | integer                                :: c_ident | 
| 138 | gezelter | 115 | integer                                :: status | 
| 139 |  |  |  | 
| 140 | chuckv | 577 | integer                                :: nAtypes,nEAMTypes,myATID | 
| 141 | gezelter | 115 | integer                                :: maxVals | 
| 142 |  |  | integer                                :: alloc_stat | 
| 143 | gezelter | 938 | integer                                :: current, j | 
| 144 | gezelter | 115 | integer,pointer                        :: Matchlist(:) => null() | 
| 145 |  |  |  | 
| 146 |  |  | status = 0 | 
| 147 |  |  |  | 
| 148 |  |  | !! Assume that atypes has already been set and get the total number of types in atypes | 
| 149 |  |  | !! Also assume that every member of atypes is a EAM model. | 
| 150 |  |  |  | 
| 151 |  |  | ! check to see if this is the first time into | 
| 152 |  |  | if (.not.associated(EAMList%EAMParams)) then | 
| 153 | chuckv | 577 | call getMatchingElementList(atypes, "is_EAM", .true., nEAMtypes, MatchList) | 
| 154 |  |  | EAMList%n_eam_types = nEAMtypes | 
| 155 |  |  | allocate(EAMList%EAMParams(nEAMTypes)) | 
| 156 |  |  | nAtypes = getSize(atypes) | 
| 157 |  |  | allocate(EAMList%atidToEAMType(nAtypes)) | 
| 158 | gezelter | 115 | end if | 
| 159 |  |  |  | 
| 160 |  |  | EAMList%currentAddition = EAMList%currentAddition + 1 | 
| 161 |  |  | current = EAMList%currentAddition | 
| 162 |  |  |  | 
| 163 | chrisfen | 579 | myATID =  getFirstMatchingElement(atypes, "c_ident", c_ident) | 
| 164 | chuckv | 577 | EAMList%atidToEAMType(myATID) = current | 
| 165 | gezelter | 507 |  | 
| 166 | chrisfen | 579 | EAMList%EAMParams(current)%eam_atype    = c_ident | 
| 167 | gezelter | 115 | EAMList%EAMParams(current)%eam_lattice  = lattice_constant | 
| 168 |  |  | EAMList%EAMParams(current)%eam_rcut     = rcut | 
| 169 |  |  |  | 
| 170 | gezelter | 938 | ! Build array of r values | 
| 171 |  |  | do j = 1, eam_nr | 
| 172 |  |  | rvals(j) = real(j-1,kind=dp) * eam_dr | 
| 173 |  |  | end do | 
| 174 |  |  |  | 
| 175 |  |  | ! Build array of rho values | 
| 176 |  |  | do j = 1, eam_nrho | 
| 177 |  |  | rhovals(j) = real(j-1,kind=dp) * eam_drho | 
| 178 |  |  | end do | 
| 179 |  |  |  | 
| 180 |  |  | ! convert from eV to kcal / mol: | 
| 181 |  |  | do j = 1, eam_nrho | 
| 182 |  |  | eam_F_rho(j) = eam_F_rho(j) * 23.06054E0_DP | 
| 183 |  |  | end do | 
| 184 |  |  |  | 
| 185 |  |  | ! precompute the pair potential and get it into kcal / mol: | 
| 186 |  |  | eam_phi_r(1) = 0.0E0_DP | 
| 187 |  |  | do j = 2, eam_nr | 
| 188 |  |  | eam_phi_r(j) = 331.999296E0_DP * (eam_Z_r(j)**2) / rvals(j) | 
| 189 |  |  | enddo | 
| 190 |  |  |  | 
| 191 | gezelter | 939 | call newSpline(EAMList%EAMParams(current)%rho, rvals, eam_rho_r, .true.) | 
| 192 |  |  | call newSpline(EAMList%EAMParams(current)%Z,   rvals, eam_Z_r, .true.) | 
| 193 |  |  | call newSpline(EAMList%EAMParams(current)%F, rhovals, eam_F_rho, .true.) | 
| 194 |  |  | call newSpline(EAMList%EAMParams(current)%phi, rvals, eam_phi_r, .true.) | 
| 195 | gezelter | 115 | end subroutine newEAMtype | 
| 196 |  |  |  | 
| 197 |  |  |  | 
| 198 | chuckv | 472 | ! kills all eam types entered and sets simulation to uninitalized | 
| 199 |  |  | subroutine destroyEAMtypes() | 
| 200 |  |  | integer :: i | 
| 201 |  |  | type(EAMType), pointer :: tempEAMType=>null() | 
| 202 | gezelter | 115 |  | 
| 203 | chuckv | 472 | do i = 1, EAMList%n_eam_types | 
| 204 |  |  | tempEAMType => eamList%EAMParams(i) | 
| 205 |  |  | call deallocate_EAMType(tempEAMType) | 
| 206 |  |  | end do | 
| 207 |  |  | if(associated( eamList%EAMParams)) deallocate( eamList%EAMParams) | 
| 208 |  |  | eamList%EAMParams => null() | 
| 209 | gezelter | 507 |  | 
| 210 | chuckv | 472 | eamList%n_eam_types = 0 | 
| 211 |  |  | eamList%currentAddition = 0 | 
| 212 |  |  | end subroutine destroyEAMtypes | 
| 213 |  |  |  | 
| 214 | chrisfen | 578 | function getEAMCut(atomID) result(cutValue) | 
| 215 |  |  | integer, intent(in) :: atomID | 
| 216 | chrisfen | 579 | integer :: eamID | 
| 217 | chrisfen | 578 | real(kind=dp) :: cutValue | 
| 218 |  |  |  | 
| 219 | chrisfen | 579 | eamID = EAMList%atidToEAMType(atomID) | 
| 220 |  |  | cutValue = EAMList%EAMParams(eamID)%eam_rcut | 
| 221 | chrisfen | 578 | end function getEAMCut | 
| 222 |  |  |  | 
| 223 | gezelter | 115 | subroutine init_EAM_FF(status) | 
| 224 |  |  | integer :: status | 
| 225 |  |  | integer :: i,j | 
| 226 |  |  | real(kind=dp) :: current_rcut_max | 
| 227 | gezelter | 938 | #ifdef IS_MPI | 
| 228 |  |  | integer :: nAtomsInRow | 
| 229 |  |  | integer :: nAtomsInCol | 
| 230 |  |  | #endif | 
| 231 | gezelter | 115 | integer :: alloc_stat | 
| 232 |  |  | integer :: number_r, number_rho | 
| 233 |  |  |  | 
| 234 |  |  | status = 0 | 
| 235 |  |  | if (EAMList%currentAddition == 0) then | 
| 236 |  |  | call handleError("init_EAM_FF","No members in EAMList") | 
| 237 |  |  | status = -1 | 
| 238 |  |  | return | 
| 239 |  |  | end if | 
| 240 | gezelter | 938 |  | 
| 241 | gezelter | 507 | !! Allocate arrays for force calculation | 
| 242 | gezelter | 938 |  | 
| 243 | gezelter | 115 | #ifdef IS_MPI | 
| 244 |  |  | nAtomsInRow = getNatomsInRow(plan_atom_row) | 
| 245 |  |  | nAtomsInCol = getNatomsInCol(plan_atom_col) | 
| 246 |  |  | #endif | 
| 247 |  |  |  | 
| 248 |  |  | if (allocated(frho)) deallocate(frho) | 
| 249 |  |  | allocate(frho(nlocal),stat=alloc_stat) | 
| 250 |  |  | if (alloc_stat /= 0) then | 
| 251 |  |  | status = -1 | 
| 252 |  |  | return | 
| 253 |  |  | end if | 
| 254 | gezelter | 938 |  | 
| 255 | gezelter | 115 | if (allocated(rho)) deallocate(rho) | 
| 256 |  |  | allocate(rho(nlocal),stat=alloc_stat) | 
| 257 |  |  | if (alloc_stat /= 0) then | 
| 258 |  |  | status = -1 | 
| 259 |  |  | return | 
| 260 |  |  | end if | 
| 261 |  |  |  | 
| 262 |  |  | if (allocated(dfrhodrho)) deallocate(dfrhodrho) | 
| 263 |  |  | allocate(dfrhodrho(nlocal),stat=alloc_stat) | 
| 264 |  |  | if (alloc_stat /= 0) then | 
| 265 |  |  | status = -1 | 
| 266 |  |  | return | 
| 267 |  |  | end if | 
| 268 |  |  |  | 
| 269 |  |  | #ifdef IS_MPI | 
| 270 |  |  |  | 
| 271 |  |  | if (allocated(rho_tmp)) deallocate(rho_tmp) | 
| 272 |  |  | allocate(rho_tmp(nlocal),stat=alloc_stat) | 
| 273 |  |  | if (alloc_stat /= 0) then | 
| 274 |  |  | status = -1 | 
| 275 |  |  | return | 
| 276 |  |  | end if | 
| 277 |  |  |  | 
| 278 |  |  | if (allocated(frho_row)) deallocate(frho_row) | 
| 279 |  |  | allocate(frho_row(nAtomsInRow),stat=alloc_stat) | 
| 280 |  |  | if (alloc_stat /= 0) then | 
| 281 |  |  | status = -1 | 
| 282 |  |  | return | 
| 283 |  |  | end if | 
| 284 |  |  | if (allocated(rho_row)) deallocate(rho_row) | 
| 285 |  |  | allocate(rho_row(nAtomsInRow),stat=alloc_stat) | 
| 286 |  |  | if (alloc_stat /= 0) then | 
| 287 |  |  | status = -1 | 
| 288 |  |  | return | 
| 289 |  |  | end if | 
| 290 |  |  | if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row) | 
| 291 |  |  | allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat) | 
| 292 |  |  | if (alloc_stat /= 0) then | 
| 293 |  |  | status = -1 | 
| 294 |  |  | return | 
| 295 |  |  | end if | 
| 296 |  |  |  | 
| 297 | gezelter | 507 | ! Now do column arrays | 
| 298 | gezelter | 115 |  | 
| 299 |  |  | if (allocated(frho_col)) deallocate(frho_col) | 
| 300 |  |  | allocate(frho_col(nAtomsInCol),stat=alloc_stat) | 
| 301 |  |  | if (alloc_stat /= 0) then | 
| 302 |  |  | status = -1 | 
| 303 |  |  | return | 
| 304 |  |  | end if | 
| 305 |  |  | if (allocated(rho_col)) deallocate(rho_col) | 
| 306 |  |  | allocate(rho_col(nAtomsInCol),stat=alloc_stat) | 
| 307 |  |  | if (alloc_stat /= 0) then | 
| 308 |  |  | status = -1 | 
| 309 |  |  | return | 
| 310 |  |  | end if | 
| 311 |  |  | if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col) | 
| 312 |  |  | allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat) | 
| 313 |  |  | if (alloc_stat /= 0) then | 
| 314 |  |  | status = -1 | 
| 315 |  |  | return | 
| 316 |  |  | end if | 
| 317 | gezelter | 507 |  | 
| 318 | gezelter | 115 | #endif | 
| 319 |  |  |  | 
| 320 | gezelter | 938 | end subroutine init_EAM_FF | 
| 321 | gezelter | 115 |  | 
| 322 | gezelter | 938 | subroutine setCutoffEAM(rcut) | 
| 323 | gezelter | 115 | real(kind=dp) :: rcut | 
| 324 |  |  | EAM_rcut = rcut | 
| 325 |  |  | end subroutine setCutoffEAM | 
| 326 |  |  |  | 
| 327 |  |  | subroutine clean_EAM() | 
| 328 | gezelter | 507 |  | 
| 329 |  |  | ! clean non-IS_MPI first | 
| 330 | gezelter | 115 | frho = 0.0_dp | 
| 331 |  |  | rho  = 0.0_dp | 
| 332 |  |  | dfrhodrho = 0.0_dp | 
| 333 | gezelter | 507 | ! clean MPI if needed | 
| 334 | gezelter | 115 | #ifdef IS_MPI | 
| 335 |  |  | frho_row = 0.0_dp | 
| 336 |  |  | frho_col = 0.0_dp | 
| 337 |  |  | rho_row  = 0.0_dp | 
| 338 |  |  | rho_col  = 0.0_dp | 
| 339 |  |  | rho_tmp  = 0.0_dp | 
| 340 |  |  | dfrhodrho_row = 0.0_dp | 
| 341 |  |  | dfrhodrho_col = 0.0_dp | 
| 342 |  |  | #endif | 
| 343 |  |  | end subroutine clean_EAM | 
| 344 |  |  |  | 
| 345 |  |  | subroutine deallocate_EAMType(thisEAMType) | 
| 346 |  |  | type (EAMtype), pointer :: thisEAMType | 
| 347 |  |  |  | 
| 348 | gezelter | 938 | call deleteSpline(thisEAMType%F) | 
| 349 |  |  | call deleteSpline(thisEAMType%rho) | 
| 350 |  |  | call deleteSpline(thisEAMType%phi) | 
| 351 |  |  | call deleteSpline(thisEAMType%Z) | 
| 352 | gezelter | 507 |  | 
| 353 | gezelter | 115 | end subroutine deallocate_EAMType | 
| 354 |  |  |  | 
| 355 | gezelter | 507 | !! Calculates rho_r | 
| 356 | gezelter | 115 | subroutine calc_eam_prepair_rho(atom1,atom2,d,r,rijsq) | 
| 357 | gezelter | 938 | integer :: atom1, atom2 | 
| 358 | gezelter | 115 | real(kind = dp), dimension(3) :: d | 
| 359 |  |  | real(kind = dp), intent(inout)               :: r | 
| 360 |  |  | real(kind = dp), intent(inout)               :: rijsq | 
| 361 |  |  | ! value of electron density rho do to atom i at atom j | 
| 362 |  |  | real(kind = dp) :: rho_i_at_j | 
| 363 |  |  | ! value of electron density rho do to atom j at atom i | 
| 364 |  |  | real(kind = dp) :: rho_j_at_i | 
| 365 |  |  | integer :: eam_err | 
| 366 | chuckv | 592 |  | 
| 367 | gezelter | 938 | integer :: atid1, atid2 ! Global atid | 
| 368 | chuckv | 592 | integer :: myid_atom1 ! EAM atid | 
| 369 |  |  | integer :: myid_atom2 | 
| 370 | gezelter | 507 |  | 
| 371 |  |  | ! check to see if we need to be cleaned at the start of a force loop | 
| 372 | gezelter | 115 |  | 
| 373 |  |  | #ifdef IS_MPI | 
| 374 | chuckv | 592 | Atid1 = Atid_row(Atom1) | 
| 375 |  |  | Atid2 = Atid_col(Atom2) | 
| 376 | gezelter | 115 | #else | 
| 377 | chuckv | 592 | Atid1 = Atid(Atom1) | 
| 378 |  |  | Atid2 = Atid(Atom2) | 
| 379 | gezelter | 115 | #endif | 
| 380 |  |  |  | 
| 381 | chuckv | 592 | Myid_atom1 = Eamlist%atidtoeamtype(Atid1) | 
| 382 |  |  | Myid_atom2 = Eamlist%atidtoeamtype(Atid2) | 
| 383 |  |  |  | 
| 384 | gezelter | 115 | if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then | 
| 385 |  |  |  | 
| 386 | chrisfen | 943 | call lookupEAMSpline(EAMList%EAMParams(myid_atom1)%rho, r, & | 
| 387 | gezelter | 938 | rho_i_at_j) | 
| 388 | gezelter | 115 |  | 
| 389 |  |  | #ifdef  IS_MPI | 
| 390 |  |  | rho_col(atom2) = rho_col(atom2) + rho_i_at_j | 
| 391 |  |  | #else | 
| 392 |  |  | rho(atom2) = rho(atom2) + rho_i_at_j | 
| 393 |  |  | #endif | 
| 394 | gezelter | 507 | endif | 
| 395 | gezelter | 115 |  | 
| 396 | gezelter | 507 | if (r.lt.EAMList%EAMParams(myid_atom2)%eam_rcut) then | 
| 397 | gezelter | 115 |  | 
| 398 | chrisfen | 943 | call lookupEAMSpline(EAMList%EAMParams(myid_atom2)%rho, r, & | 
| 399 | gezelter | 938 | rho_j_at_i) | 
| 400 | gezelter | 507 |  | 
| 401 | gezelter | 115 | #ifdef  IS_MPI | 
| 402 | gezelter | 507 | rho_row(atom1) = rho_row(atom1) + rho_j_at_i | 
| 403 | gezelter | 115 | #else | 
| 404 | gezelter | 507 | rho(atom1) = rho(atom1) + rho_j_at_i | 
| 405 | gezelter | 115 | #endif | 
| 406 | gezelter | 507 | endif | 
| 407 | gezelter | 115 | end subroutine calc_eam_prepair_rho | 
| 408 |  |  |  | 
| 409 |  |  |  | 
| 410 |  |  | !! Calculate the functional F(rho) for all local atoms | 
| 411 | gezelter | 938 | subroutine calc_eam_preforce_Frho(nlocal, pot) | 
| 412 | gezelter | 115 | integer :: nlocal | 
| 413 |  |  | real(kind=dp) :: pot | 
| 414 | gezelter | 938 | integer :: i, j | 
| 415 | gezelter | 115 | integer :: atom | 
| 416 | gezelter | 938 | real(kind=dp) :: U,U1 | 
| 417 | gezelter | 115 | integer :: atype1 | 
| 418 | gezelter | 938 | integer :: me, atid1 | 
| 419 | gezelter | 115 |  | 
| 420 |  |  | cleanme = .true. | 
| 421 | gezelter | 938 | !! Scatter the electron density from  pre-pair calculation back to | 
| 422 |  |  | !! local atoms | 
| 423 | gezelter | 115 | #ifdef IS_MPI | 
| 424 |  |  | call scatter(rho_row,rho,plan_atom_row,eam_err) | 
| 425 |  |  | if (eam_err /= 0 ) then | 
| 426 | gezelter | 507 | write(errMsg,*) " Error scattering rho_row into rho" | 
| 427 |  |  | call handleError(RoutineName,errMesg) | 
| 428 |  |  | endif | 
| 429 | gezelter | 115 | call scatter(rho_col,rho_tmp,plan_atom_col,eam_err) | 
| 430 |  |  | if (eam_err /= 0 ) then | 
| 431 | gezelter | 507 | write(errMsg,*) " Error scattering rho_col into rho" | 
| 432 |  |  | call handleError(RoutineName,errMesg) | 
| 433 |  |  | endif | 
| 434 | gezelter | 115 |  | 
| 435 | gezelter | 507 | rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal) | 
| 436 | gezelter | 115 | #endif | 
| 437 |  |  |  | 
| 438 | gezelter | 507 | !! Calculate F(rho) and derivative | 
| 439 | gezelter | 115 | do atom = 1, nlocal | 
| 440 | chuckv | 592 | atid1 = atid(atom) | 
| 441 |  |  | me = eamList%atidToEAMtype(atid1) | 
| 442 | gezelter | 115 |  | 
| 443 | chrisfen | 943 | call lookupEAMSpline1d(EAMList%EAMParams(me)%F, rho(atom), & | 
| 444 | gezelter | 938 | u, u1) | 
| 445 |  |  |  | 
| 446 | gezelter | 115 | frho(atom) = u | 
| 447 |  |  | dfrhodrho(atom) = u1 | 
| 448 |  |  | pot = pot + u | 
| 449 | gezelter | 939 |  | 
| 450 | gezelter | 115 | enddo | 
| 451 |  |  |  | 
| 452 |  |  | #ifdef IS_MPI | 
| 453 |  |  | !! communicate f(rho) and derivatives back into row and column arrays | 
| 454 |  |  | call gather(frho,frho_row,plan_atom_row, eam_err) | 
| 455 |  |  | if (eam_err /=  0) then | 
| 456 |  |  | call handleError("cal_eam_forces()","MPI gather frho_row failure") | 
| 457 |  |  | endif | 
| 458 |  |  | call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, eam_err) | 
| 459 |  |  | if (eam_err /=  0) then | 
| 460 |  |  | call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure") | 
| 461 |  |  | endif | 
| 462 |  |  | call gather(frho,frho_col,plan_atom_col, eam_err) | 
| 463 |  |  | if (eam_err /=  0) then | 
| 464 |  |  | call handleError("cal_eam_forces()","MPI gather frho_col failure") | 
| 465 |  |  | endif | 
| 466 |  |  | call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, eam_err) | 
| 467 |  |  | if (eam_err /=  0) then | 
| 468 |  |  | call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure") | 
| 469 |  |  | endif | 
| 470 |  |  | #endif | 
| 471 |  |  |  | 
| 472 | gezelter | 507 |  | 
| 473 | gezelter | 115 | end subroutine calc_eam_preforce_Frho | 
| 474 | gezelter | 938 |  | 
| 475 | gezelter | 115 | !! Does EAM pairwise Force calculation. | 
| 476 |  |  | subroutine do_eam_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, & | 
| 477 |  |  | pot, f, do_pot) | 
| 478 |  |  | !Arguments | 
| 479 |  |  | integer, intent(in) ::  atom1, atom2 | 
| 480 |  |  | real( kind = dp ), intent(in) :: rij, r2 | 
| 481 |  |  | real( kind = dp ) :: pot, sw, vpair | 
| 482 |  |  | real( kind = dp ), dimension(3,nLocal) :: f | 
| 483 |  |  | real( kind = dp ), intent(in), dimension(3) :: d | 
| 484 |  |  | real( kind = dp ), intent(inout), dimension(3) :: fpair | 
| 485 |  |  |  | 
| 486 |  |  | logical, intent(in) :: do_pot | 
| 487 | gezelter | 507 |  | 
| 488 | gezelter | 938 | real( kind = dp ) :: drdx, drdy, drdz | 
| 489 |  |  | real( kind = dp ) :: phab, pha, dvpdr | 
| 490 |  |  | real( kind = dp ) :: rha, drha, dpha | 
| 491 |  |  | real( kind = dp ) :: rhb, drhb, dphb | 
| 492 | gezelter | 115 | real( kind = dp ) :: dudr | 
| 493 | gezelter | 938 | real( kind = dp ) :: rci, rcj | 
| 494 |  |  | real( kind = dp ) :: drhoidr, drhojdr | 
| 495 |  |  | real( kind = dp ) :: Fx, Fy, Fz | 
| 496 |  |  | real( kind = dp ) :: r, phb | 
| 497 | gezelter | 115 |  | 
| 498 | gezelter | 938 | integer :: id1, id2 | 
| 499 | gezelter | 115 | integer  :: mytype_atom1 | 
| 500 |  |  | integer  :: mytype_atom2 | 
| 501 | gezelter | 938 | integer  :: atid1, atid2 | 
| 502 | gezelter | 115 |  | 
| 503 |  |  | phab = 0.0E0_DP | 
| 504 |  |  | dvpdr = 0.0E0_DP | 
| 505 |  |  |  | 
| 506 |  |  | if (rij .lt. EAM_rcut) then | 
| 507 |  |  |  | 
| 508 |  |  | #ifdef IS_MPI | 
| 509 | chuckv | 592 | atid1 = atid_row(atom1) | 
| 510 |  |  | atid2 = atid_col(atom2) | 
| 511 | gezelter | 115 | #else | 
| 512 | chuckv | 592 | atid1 = atid(atom1) | 
| 513 |  |  | atid2 = atid(atom2) | 
| 514 | gezelter | 115 | #endif | 
| 515 | chuckv | 592 |  | 
| 516 |  |  | mytype_atom1 = EAMList%atidToEAMType(atid1) | 
| 517 |  |  | mytype_atom2 = EAMList%atidTOEAMType(atid2) | 
| 518 |  |  |  | 
| 519 |  |  |  | 
| 520 | gezelter | 115 | ! get cutoff for atom 1 | 
| 521 |  |  | rci = EAMList%EAMParams(mytype_atom1)%eam_rcut | 
| 522 |  |  | ! get type specific cutoff for atom 2 | 
| 523 |  |  | rcj = EAMList%EAMParams(mytype_atom2)%eam_rcut | 
| 524 | gezelter | 507 |  | 
| 525 | gezelter | 115 | drdx = d(1)/rij | 
| 526 |  |  | drdy = d(2)/rij | 
| 527 |  |  | drdz = d(3)/rij | 
| 528 | gezelter | 507 |  | 
| 529 | gezelter | 115 | if (rij.lt.rci) then | 
| 530 | gezelter | 938 |  | 
| 531 |  |  | ! Calculate rho and drho for atom1 | 
| 532 |  |  |  | 
| 533 | chrisfen | 943 | call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%rho, & | 
| 534 | gezelter | 938 | rij, rha, drha) | 
| 535 |  |  |  | 
| 536 |  |  | ! Calculate Phi(r) for atom1. | 
| 537 |  |  |  | 
| 538 | chrisfen | 943 | call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom1)%phi, & | 
| 539 | gezelter | 938 | rij, pha, dpha) | 
| 540 |  |  |  | 
| 541 | gezelter | 115 | endif | 
| 542 |  |  |  | 
| 543 |  |  | if (rij.lt.rcj) then | 
| 544 | gezelter | 507 |  | 
| 545 | gezelter | 938 | ! Calculate rho and drho for atom2 | 
| 546 |  |  |  | 
| 547 | chrisfen | 943 | call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%rho, & | 
| 548 | gezelter | 938 | rij, rhb, drhb) | 
| 549 |  |  |  | 
| 550 |  |  | ! Calculate Phi(r) for atom2. | 
| 551 |  |  |  | 
| 552 | chrisfen | 943 | call lookupEAMSpline1d(EAMList%EAMParams(mytype_atom2)%phi, & | 
| 553 | gezelter | 938 | rij, phb, dphb) | 
| 554 |  |  |  | 
| 555 | gezelter | 115 | endif | 
| 556 |  |  |  | 
| 557 |  |  | if (rij.lt.rci) then | 
| 558 |  |  | phab = phab + 0.5E0_DP*(rhb/rha)*pha | 
| 559 |  |  | dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + & | 
| 560 |  |  | pha*((drhb/rha) - (rhb*drha/rha/rha))) | 
| 561 |  |  | endif | 
| 562 | gezelter | 507 |  | 
| 563 | gezelter | 115 | if (rij.lt.rcj) then | 
| 564 |  |  | phab = phab + 0.5E0_DP*(rha/rhb)*phb | 
| 565 |  |  | dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + & | 
| 566 |  |  | phb*((drha/rhb) - (rha*drhb/rhb/rhb))) | 
| 567 |  |  | endif | 
| 568 | gezelter | 507 |  | 
| 569 | gezelter | 115 | drhoidr = drha | 
| 570 |  |  | drhojdr = drhb | 
| 571 |  |  |  | 
| 572 |  |  | #ifdef IS_MPI | 
| 573 |  |  | dudr = drhojdr*dfrhodrho_row(atom1)+drhoidr*dfrhodrho_col(atom2) & | 
| 574 |  |  | + dvpdr | 
| 575 |  |  |  | 
| 576 |  |  | #else | 
| 577 |  |  | dudr = drhojdr*dfrhodrho(atom1)+drhoidr*dfrhodrho(atom2) & | 
| 578 |  |  | + dvpdr | 
| 579 |  |  | #endif | 
| 580 |  |  |  | 
| 581 |  |  | fx = dudr * drdx | 
| 582 |  |  | fy = dudr * drdy | 
| 583 |  |  | fz = dudr * drdz | 
| 584 |  |  |  | 
| 585 |  |  |  | 
| 586 |  |  | #ifdef IS_MPI | 
| 587 |  |  | if (do_pot) then | 
| 588 | gezelter | 662 | pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + phab*0.5 | 
| 589 |  |  | pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + phab*0.5 | 
| 590 | gezelter | 115 | end if | 
| 591 |  |  |  | 
| 592 |  |  | f_Row(1,atom1) = f_Row(1,atom1) + fx | 
| 593 |  |  | f_Row(2,atom1) = f_Row(2,atom1) + fy | 
| 594 |  |  | f_Row(3,atom1) = f_Row(3,atom1) + fz | 
| 595 | gezelter | 507 |  | 
| 596 | gezelter | 115 | f_Col(1,atom2) = f_Col(1,atom2) - fx | 
| 597 |  |  | f_Col(2,atom2) = f_Col(2,atom2) - fy | 
| 598 |  |  | f_Col(3,atom2) = f_Col(3,atom2) - fz | 
| 599 |  |  | #else | 
| 600 |  |  |  | 
| 601 |  |  | if(do_pot) then | 
| 602 |  |  | pot = pot + phab | 
| 603 |  |  | end if | 
| 604 |  |  |  | 
| 605 |  |  | f(1,atom1) = f(1,atom1) + fx | 
| 606 |  |  | f(2,atom1) = f(2,atom1) + fy | 
| 607 |  |  | f(3,atom1) = f(3,atom1) + fz | 
| 608 | gezelter | 507 |  | 
| 609 | gezelter | 115 | f(1,atom2) = f(1,atom2) - fx | 
| 610 |  |  | f(2,atom2) = f(2,atom2) - fy | 
| 611 |  |  | f(3,atom2) = f(3,atom2) - fz | 
| 612 |  |  | #endif | 
| 613 | gezelter | 507 |  | 
| 614 | gezelter | 115 | vpair = vpair + phab | 
| 615 | gezelter | 938 |  | 
| 616 | gezelter | 115 | #ifdef IS_MPI | 
| 617 |  |  | id1 = AtomRowToGlobal(atom1) | 
| 618 |  |  | id2 = AtomColToGlobal(atom2) | 
| 619 |  |  | #else | 
| 620 |  |  | id1 = atom1 | 
| 621 |  |  | id2 = atom2 | 
| 622 |  |  | #endif | 
| 623 | gezelter | 507 |  | 
| 624 | gezelter | 115 | if (molMembershipList(id1) .ne. molMembershipList(id2)) then | 
| 625 | gezelter | 507 |  | 
| 626 | gezelter | 115 | fpair(1) = fpair(1) + fx | 
| 627 |  |  | fpair(2) = fpair(2) + fy | 
| 628 |  |  | fpair(3) = fpair(3) + fz | 
| 629 | gezelter | 507 |  | 
| 630 | gezelter | 115 | endif | 
| 631 | gezelter | 507 | endif | 
| 632 | gezelter | 115 | end subroutine do_eam_pair | 
| 633 |  |  |  | 
| 634 | chrisfen | 943 | subroutine lookupEAMSpline(cs, xval, yval) | 
| 635 |  |  |  | 
| 636 |  |  | implicit none | 
| 637 |  |  |  | 
| 638 |  |  | type (cubicSpline), intent(in) :: cs | 
| 639 |  |  | real( kind = DP ), intent(in)  :: xval | 
| 640 |  |  | real( kind = DP ), intent(out) :: yval | 
| 641 |  |  | real( kind = DP ) ::  dx | 
| 642 |  |  | integer :: i, j | 
| 643 |  |  | ! | 
| 644 |  |  | !  Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains | 
| 645 |  |  | !  or is nearest to xval. | 
| 646 |  |  |  | 
| 647 |  |  | j = MAX(1, MIN(cs%n-1, idint((xval-cs%x(1)) * cs%dx_i) + 1)) | 
| 648 |  |  |  | 
| 649 |  |  | dx = xval - cs%x(j) | 
| 650 |  |  | yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j))) | 
| 651 |  |  |  | 
| 652 |  |  | return | 
| 653 |  |  | end subroutine lookupEAMSpline | 
| 654 |  |  |  | 
| 655 |  |  | subroutine lookupEAMSpline1d(cs, xval, yval, dydx) | 
| 656 |  |  |  | 
| 657 |  |  | implicit none | 
| 658 |  |  |  | 
| 659 |  |  | type (cubicSpline), intent(in) :: cs | 
| 660 |  |  | real( kind = DP ), intent(in)  :: xval | 
| 661 |  |  | real( kind = DP ), intent(out) :: yval, dydx | 
| 662 |  |  | real( kind = DP ) :: dx | 
| 663 |  |  | integer :: i, j | 
| 664 |  |  |  | 
| 665 |  |  | !  Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains | 
| 666 |  |  | !  or is nearest to xval. | 
| 667 |  |  |  | 
| 668 |  |  |  | 
| 669 |  |  | j = MAX(1, MIN(cs%n-1, idint((xval-cs%x(1)) * cs%dx_i) + 1)) | 
| 670 |  |  |  | 
| 671 |  |  | dx = xval - cs%x(j) | 
| 672 |  |  | yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j))) | 
| 673 |  |  |  | 
| 674 |  |  | dydx = cs%b(j) + dx*(2.0d0 * cs%c(j) + 3.0d0 * dx * cs%d(j)) | 
| 675 |  |  |  | 
| 676 |  |  | return | 
| 677 |  |  | end subroutine lookupEAMSpline1d | 
| 678 |  |  |  | 
| 679 | gezelter | 115 | end module eam |