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