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