ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/UseTheForce/DarkSide/eam.F90
Revision: 2276
Committed: Fri Aug 26 20:34:24 2005 UTC (18 years, 10 months ago) by chuckv
File size: 32394 byte(s)
Log Message:
Added eamType map to atid map.

File Contents

# Content
1 !!
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 module eam
43 use simulation
44 use force_globals
45 use status
46 use atype_module
47 use Vector_class
48 #ifdef IS_MPI
49 use mpiSimulation
50 #endif
51 implicit none
52 PRIVATE
53
54 INTEGER, PARAMETER :: DP = selected_real_kind(15)
55
56 logical, save :: EAM_FF_initialized = .false.
57 integer, save :: EAM_Mixing_Policy
58 real(kind = dp), save :: EAM_rcut
59 logical, save :: haveRcut = .false.
60
61 character(len = statusMsgSize) :: errMesg
62 integer :: eam_err
63
64 character(len = 200) :: errMsg
65 character(len=*), parameter :: RoutineName = "EAM MODULE"
66 !! Logical that determines if eam arrays should be zeroed
67 logical :: cleanme = .true.
68 logical :: nmflag = .false.
69
70
71 type, private :: EAMtype
72 integer :: eam_atype
73 real( kind = DP ) :: eam_dr
74 integer :: eam_nr
75 integer :: eam_nrho
76 real( kind = DP ) :: eam_lattice
77 real( kind = DP ) :: eam_drho
78 real( kind = DP ) :: eam_rcut
79 integer :: eam_atype_map
80
81 real( kind = DP ), pointer, dimension(:) :: eam_rvals => null()
82 real( kind = DP ), pointer, dimension(:) :: eam_rhovals => null()
83 real( kind = DP ), pointer, dimension(:) :: eam_F_rho => null()
84 real( kind = DP ), pointer, dimension(:) :: eam_Z_r => null()
85 real( kind = DP ), pointer, dimension(:) :: eam_rho_r => null()
86 real( kind = DP ), pointer, dimension(:) :: eam_phi_r => null()
87 real( kind = DP ), pointer, dimension(:) :: eam_F_rho_pp => null()
88 real( kind = DP ), pointer, dimension(:) :: eam_Z_r_pp => null()
89 real( kind = DP ), pointer, dimension(:) :: eam_rho_r_pp => null()
90 real( kind = DP ), pointer, dimension(:) :: eam_phi_r_pp => null()
91 end type EAMtype
92
93
94 !! Arrays for derivatives used in force calculation
95 real( kind = dp), dimension(:), allocatable :: frho
96 real( kind = dp), dimension(:), allocatable :: rho
97
98 real( kind = dp), dimension(:), allocatable :: dfrhodrho
99 real( kind = dp), dimension(:), allocatable :: d2frhodrhodrho
100
101
102 !! Arrays for MPI storage
103 #ifdef IS_MPI
104 real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_col
105 real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_row
106 real( kind = dp),save, dimension(:), allocatable :: frho_row
107 real( kind = dp),save, dimension(:), allocatable :: frho_col
108 real( kind = dp),save, dimension(:), allocatable :: rho_row
109 real( kind = dp),save, dimension(:), allocatable :: rho_col
110 real( kind = dp),save, dimension(:), allocatable :: rho_tmp
111 real( kind = dp),save, dimension(:), allocatable :: d2frhodrhodrho_col
112 real( kind = dp),save, dimension(:), allocatable :: d2frhodrhodrho_row
113 #endif
114
115 type, private :: EAMTypeList
116 integer :: n_eam_types = 0
117 integer :: currentAddition = 0
118
119 type (EAMtype), pointer :: EAMParams(:) => null()
120 integer, pointer :: atidToEAMType(:) => null()
121 end type EAMTypeList
122
123
124 type (eamTypeList), save :: EAMList
125
126 !! standard eam stuff
127
128
129 public :: init_EAM_FF
130 public :: setCutoffEAM
131 public :: do_eam_pair
132 public :: newEAMtype
133 public :: calc_eam_prepair_rho
134 public :: calc_eam_preforce_Frho
135 public :: clean_EAM
136 public :: destroyEAMTypes
137
138 contains
139
140
141 subroutine newEAMtype(lattice_constant,eam_nrho,eam_drho,eam_nr,&
142 eam_dr,rcut,eam_Z_r,eam_rho_r,eam_F_rho,&
143 eam_ident,status)
144 real (kind = dp ) :: lattice_constant
145 integer :: eam_nrho
146 real (kind = dp ) :: eam_drho
147 integer :: eam_nr
148 real (kind = dp ) :: eam_dr
149 real (kind = dp ) :: rcut
150 real (kind = dp ), dimension(eam_nr) :: eam_Z_r
151 real (kind = dp ), dimension(eam_nr) :: eam_rho_r
152 real (kind = dp ), dimension(eam_nrho) :: eam_F_rho
153 integer :: eam_ident
154 integer :: status
155
156 integer :: nAtypes,nEAMTypes,myATID
157 integer :: maxVals
158 integer :: alloc_stat
159 integer :: current
160 integer,pointer :: Matchlist(:) => null()
161
162 status = 0
163
164
165 !! Assume that atypes has already been set and get the total number of types in atypes
166 !! Also assume that every member of atypes is a EAM model.
167
168
169 ! check to see if this is the first time into
170 if (.not.associated(EAMList%EAMParams)) then
171 call getMatchingElementList(atypes, "is_EAM", .true., nEAMtypes, MatchList)
172 EAMList%n_eam_types = nEAMtypes
173 allocate(EAMList%EAMParams(nEAMTypes))
174 nAtypes = getSize(atypes)
175 allocate(EAMList%atidToEAMType(nAtypes))
176 end if
177
178 EAMList%currentAddition = EAMList%currentAddition + 1
179 current = EAMList%currentAddition
180
181 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
182 EAMList%atidToEAMType(myATID) = current
183
184 call allocate_EAMType(eam_nrho,eam_nr,EAMList%EAMParams(current),stat=alloc_stat)
185 if (alloc_stat /= 0) then
186 status = -1
187 return
188 end if
189
190
191 EAMList%EAMParams(current)%eam_atype = eam_ident
192 EAMList%EAMParams(current)%eam_lattice = lattice_constant
193 EAMList%EAMParams(current)%eam_nrho = eam_nrho
194 EAMList%EAMParams(current)%eam_drho = eam_drho
195 EAMList%EAMParams(current)%eam_nr = eam_nr
196 EAMList%EAMParams(current)%eam_dr = eam_dr
197 EAMList%EAMParams(current)%eam_rcut = rcut
198 EAMList%EAMParams(current)%eam_Z_r = eam_Z_r
199 EAMList%EAMParams(current)%eam_rho_r = eam_rho_r
200 EAMList%EAMParams(current)%eam_F_rho = eam_F_rho
201
202 end subroutine newEAMtype
203
204
205 ! kills all eam types entered and sets simulation to uninitalized
206 subroutine destroyEAMtypes()
207 integer :: i
208 type(EAMType), pointer :: tempEAMType=>null()
209
210 do i = 1, EAMList%n_eam_types
211 tempEAMType => eamList%EAMParams(i)
212 call deallocate_EAMType(tempEAMType)
213 end do
214 if(associated( eamList%EAMParams)) deallocate( eamList%EAMParams)
215 eamList%EAMParams => null()
216
217 eamList%n_eam_types = 0
218 eamList%currentAddition = 0
219
220 end subroutine destroyEAMtypes
221
222 subroutine init_EAM_FF(status)
223 integer :: status
224 integer :: i,j
225 real(kind=dp) :: current_rcut_max
226 integer :: alloc_stat
227 integer :: number_r, number_rho
228
229
230 status = 0
231 if (EAMList%currentAddition == 0) then
232 call handleError("init_EAM_FF","No members in EAMList")
233 status = -1
234 return
235 end if
236
237
238 do i = 1, EAMList%currentAddition
239
240 ! Build array of r values
241
242 do j = 1,EAMList%EAMParams(i)%eam_nr
243 EAMList%EAMParams(i)%eam_rvals(j) = &
244 real(j-1,kind=dp)* &
245 EAMList%EAMParams(i)%eam_dr
246 end do
247 ! Build array of rho values
248 do j = 1,EAMList%EAMParams(i)%eam_nrho
249 EAMList%EAMParams(i)%eam_rhovals(j) = &
250 real(j-1,kind=dp)* &
251 EAMList%EAMParams(i)%eam_drho
252 end do
253 ! convert from eV to kcal / mol:
254 EAMList%EAMParams(i)%eam_F_rho = EAMList%EAMParams(i)%eam_F_rho * 23.06054E0_DP
255
256 ! precompute the pair potential and get it into kcal / mol:
257 EAMList%EAMParams(i)%eam_phi_r(1) = 0.0E0_DP
258 do j = 2, EAMList%EAMParams(i)%eam_nr
259 EAMList%EAMParams(i)%eam_phi_r(j) = (EAMList%EAMParams(i)%eam_Z_r(j)**2)/EAMList%EAMParams(i)%eam_rvals(j)
260 EAMList%EAMParams(i)%eam_phi_r(j) = EAMList%EAMParams(i)%eam_phi_r(j)*331.999296E0_DP
261 enddo
262 end do
263
264
265 do i = 1, EAMList%currentAddition
266 number_r = EAMList%EAMParams(i)%eam_nr
267 number_rho = EAMList%EAMParams(i)%eam_nrho
268
269 call eam_spline(number_r, EAMList%EAMParams(i)%eam_rvals, &
270 EAMList%EAMParams(i)%eam_rho_r, &
271 EAMList%EAMParams(i)%eam_rho_r_pp, &
272 0.0E0_DP, 0.0E0_DP, 'N')
273 call eam_spline(number_r, EAMList%EAMParams(i)%eam_rvals, &
274 EAMList%EAMParams(i)%eam_Z_r, &
275 EAMList%EAMParams(i)%eam_Z_r_pp, &
276 0.0E0_DP, 0.0E0_DP, 'N')
277 call eam_spline(number_rho, EAMList%EAMParams(i)%eam_rhovals, &
278 EAMList%EAMParams(i)%eam_F_rho, &
279 EAMList%EAMParams(i)%eam_F_rho_pp, &
280 0.0E0_DP, 0.0E0_DP, 'N')
281 call eam_spline(number_r, EAMList%EAMParams(i)%eam_rvals, &
282 EAMList%EAMParams(i)%eam_phi_r, &
283 EAMList%EAMParams(i)%eam_phi_r_pp, &
284 0.0E0_DP, 0.0E0_DP, 'N')
285 enddo
286
287 ! current_rcut_max = EAMList%EAMParams(1)%eam_rcut
288 !! find the smallest rcut for any eam atype
289 ! do i = 2, EAMList%currentAddition
290 ! current_rcut_max =max(current_rcut_max,EAMList%EAMParams(i)%eam_rcut)
291 ! end do
292
293 ! EAM_rcut = current_rcut_max
294 ! EAM_rcut_orig = current_rcut_max
295 ! do i = 1, EAMList%currentAddition
296 ! EAMList%EAMParam(i)s%eam_atype_map(eam_atype(i)) = i
297 ! end do
298 !! Allocate arrays for force calculation
299
300 call allocateEAM(alloc_stat)
301 if (alloc_stat /= 0 ) then
302 write(*,*) "allocateEAM failed"
303 status = -1
304 return
305 endif
306
307 end subroutine init_EAM_FF
308
309 !! routine checks to see if array is allocated, deallocates array if allocated
310 !! and then creates the array to the required size
311 subroutine allocateEAM(status)
312 integer, intent(out) :: status
313
314 #ifdef IS_MPI
315 integer :: nAtomsInRow
316 integer :: nAtomsInCol
317 #endif
318 integer :: alloc_stat
319
320
321 status = 0
322 #ifdef IS_MPI
323 nAtomsInRow = getNatomsInRow(plan_atom_row)
324 nAtomsInCol = getNatomsInCol(plan_atom_col)
325 #endif
326
327 if (allocated(frho)) deallocate(frho)
328 allocate(frho(nlocal),stat=alloc_stat)
329 if (alloc_stat /= 0) then
330 status = -1
331 return
332 end if
333 if (allocated(rho)) deallocate(rho)
334 allocate(rho(nlocal),stat=alloc_stat)
335 if (alloc_stat /= 0) then
336 status = -1
337 return
338 end if
339
340 if (allocated(dfrhodrho)) deallocate(dfrhodrho)
341 allocate(dfrhodrho(nlocal),stat=alloc_stat)
342 if (alloc_stat /= 0) then
343 status = -1
344 return
345 end if
346
347 if (allocated(d2frhodrhodrho)) deallocate(d2frhodrhodrho)
348 allocate(d2frhodrhodrho(nlocal),stat=alloc_stat)
349 if (alloc_stat /= 0) then
350 status = -1
351 return
352 end if
353
354 #ifdef IS_MPI
355
356 if (allocated(rho_tmp)) deallocate(rho_tmp)
357 allocate(rho_tmp(nlocal),stat=alloc_stat)
358 if (alloc_stat /= 0) then
359 status = -1
360 return
361 end if
362
363
364 if (allocated(frho_row)) deallocate(frho_row)
365 allocate(frho_row(nAtomsInRow),stat=alloc_stat)
366 if (alloc_stat /= 0) then
367 status = -1
368 return
369 end if
370 if (allocated(rho_row)) deallocate(rho_row)
371 allocate(rho_row(nAtomsInRow),stat=alloc_stat)
372 if (alloc_stat /= 0) then
373 status = -1
374 return
375 end if
376 if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row)
377 allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat)
378 if (alloc_stat /= 0) then
379 status = -1
380 return
381 end if
382 if (allocated(d2frhodrhodrho_row)) deallocate(d2frhodrhodrho_row)
383 allocate(d2frhodrhodrho_row(nAtomsInRow),stat=alloc_stat)
384 if (alloc_stat /= 0) then
385 status = -1
386 return
387 end if
388
389
390 ! Now do column arrays
391
392 if (allocated(frho_col)) deallocate(frho_col)
393 allocate(frho_col(nAtomsInCol),stat=alloc_stat)
394 if (alloc_stat /= 0) then
395 status = -1
396 return
397 end if
398 if (allocated(rho_col)) deallocate(rho_col)
399 allocate(rho_col(nAtomsInCol),stat=alloc_stat)
400 if (alloc_stat /= 0) then
401 status = -1
402 return
403 end if
404 if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col)
405 allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat)
406 if (alloc_stat /= 0) then
407 status = -1
408 return
409 end if
410 if (allocated(d2frhodrhodrho_col)) deallocate(d2frhodrhodrho_col)
411 allocate(d2frhodrhodrho_col(nAtomsInCol),stat=alloc_stat)
412 if (alloc_stat /= 0) then
413 status = -1
414 return
415 end if
416
417 #endif
418
419 end subroutine allocateEAM
420
421 !! C sets rcut to be the largest cutoff of any atype
422 !! present in this simulation. Doesn't include all atypes
423 !! sim knows about, just those in the simulation.
424 subroutine setCutoffEAM(rcut, status)
425 real(kind=dp) :: rcut
426 integer :: status
427 status = 0
428
429 EAM_rcut = rcut
430
431 end subroutine setCutoffEAM
432
433
434
435 subroutine clean_EAM()
436
437 ! clean non-IS_MPI first
438 frho = 0.0_dp
439 rho = 0.0_dp
440 dfrhodrho = 0.0_dp
441 ! clean MPI if needed
442 #ifdef IS_MPI
443 frho_row = 0.0_dp
444 frho_col = 0.0_dp
445 rho_row = 0.0_dp
446 rho_col = 0.0_dp
447 rho_tmp = 0.0_dp
448 dfrhodrho_row = 0.0_dp
449 dfrhodrho_col = 0.0_dp
450 #endif
451 end subroutine clean_EAM
452
453
454
455 subroutine allocate_EAMType(eam_n_rho,eam_n_r,thisEAMType,stat)
456 integer, intent(in) :: eam_n_rho
457 integer, intent(in) :: eam_n_r
458 type (EAMType) :: thisEAMType
459 integer, optional :: stat
460 integer :: alloc_stat
461
462
463
464 if (present(stat)) stat = 0
465
466 allocate(thisEAMType%eam_rvals(eam_n_r),stat=alloc_stat)
467 if (alloc_stat /= 0 ) then
468 if (present(stat)) stat = -1
469 return
470 end if
471 allocate(thisEAMType%eam_rhovals(eam_n_rho),stat=alloc_stat)
472 if (alloc_stat /= 0 ) then
473 if (present(stat)) stat = -1
474 return
475 end if
476 allocate(thisEAMType%eam_F_rho(eam_n_rho),stat=alloc_stat)
477 if (alloc_stat /= 0 ) then
478 if (present(stat)) stat = -1
479 return
480 end if
481 allocate(thisEAMType%eam_Z_r(eam_n_r),stat=alloc_stat)
482 if (alloc_stat /= 0 ) then
483 if (present(stat)) stat = -1
484 return
485 end if
486 allocate(thisEAMType%eam_rho_r(eam_n_r),stat=alloc_stat)
487 if (alloc_stat /= 0 ) then
488 if (present(stat)) stat = -1
489 return
490 end if
491 allocate(thisEAMType%eam_phi_r(eam_n_r),stat=alloc_stat)
492 if (alloc_stat /= 0 ) then
493 if (present(stat)) stat = -1
494 return
495 end if
496 allocate(thisEAMType%eam_F_rho_pp(eam_n_rho),stat=alloc_stat)
497 if (alloc_stat /= 0 ) then
498 if (present(stat)) stat = -1
499 return
500 end if
501 allocate(thisEAMType%eam_Z_r_pp(eam_n_r),stat=alloc_stat)
502 if (alloc_stat /= 0 ) then
503 if (present(stat)) stat = -1
504 return
505 end if
506 allocate(thisEAMType%eam_rho_r_pp(eam_n_r),stat=alloc_stat)
507 if (alloc_stat /= 0 ) then
508 if (present(stat)) stat = -1
509 return
510 end if
511 allocate(thisEAMType%eam_phi_r_pp(eam_n_r),stat=alloc_stat)
512 if (alloc_stat /= 0 ) then
513 if (present(stat)) stat = -1
514 return
515 end if
516
517
518 end subroutine allocate_EAMType
519
520
521 subroutine deallocate_EAMType(thisEAMType)
522 type (EAMtype), pointer :: thisEAMType
523
524 ! free Arrays in reverse order of allocation...
525 if(associated(thisEAMType%eam_phi_r_pp)) deallocate(thisEAMType%eam_phi_r_pp)
526 if(associated(thisEAMType%eam_rho_r_pp)) deallocate(thisEAMType%eam_rho_r_pp)
527 if(associated(thisEAMType%eam_Z_r_pp)) deallocate(thisEAMType%eam_Z_r_pp)
528 if(associated(thisEAMType%eam_F_rho_pp)) deallocate(thisEAMType%eam_F_rho_pp)
529 if(associated(thisEAMType%eam_phi_r)) deallocate(thisEAMType%eam_phi_r)
530 if(associated(thisEAMType%eam_rho_r)) deallocate(thisEAMType%eam_rho_r)
531 if(associated(thisEAMType%eam_Z_r)) deallocate(thisEAMType%eam_Z_r)
532 if(associated(thisEAMType%eam_F_rho)) deallocate(thisEAMType%eam_F_rho)
533 if(associated(thisEAMType%eam_rhovals)) deallocate(thisEAMType%eam_rhovals)
534 if(associated(thisEAMType%eam_rvals)) deallocate(thisEAMType%eam_rvals)
535
536 end subroutine deallocate_EAMType
537
538 !! Calculates rho_r
539 subroutine calc_eam_prepair_rho(atom1,atom2,d,r,rijsq)
540 integer :: atom1,atom2
541 real(kind = dp), dimension(3) :: d
542 real(kind = dp), intent(inout) :: r
543 real(kind = dp), intent(inout) :: rijsq
544 ! value of electron density rho do to atom i at atom j
545 real(kind = dp) :: rho_i_at_j
546 ! value of electron density rho do to atom j at atom i
547 real(kind = dp) :: rho_j_at_i
548
549 ! we don't use the derivatives, dummy variables
550 real( kind = dp) :: drho,d2rho
551 integer :: eam_err
552
553 integer :: myid_atom1
554 integer :: myid_atom2
555
556 ! check to see if we need to be cleaned at the start of a force loop
557
558
559
560
561 #ifdef IS_MPI
562 myid_atom1 = atid_Row(atom1)
563 myid_atom2 = atid_Col(atom2)
564 #else
565 myid_atom1 = atid(atom1)
566 myid_atom2 = atid(atom2)
567 #endif
568
569 if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
570
571
572
573 call eam_splint(EAMList%EAMParams(myid_atom1)%eam_nr, &
574 EAMList%EAMParams(myid_atom1)%eam_rvals, &
575 EAMList%EAMParams(myid_atom1)%eam_rho_r, &
576 EAMList%EAMParams(myid_atom1)%eam_rho_r_pp, &
577 r, rho_i_at_j,drho,d2rho)
578
579
580
581 #ifdef IS_MPI
582 rho_col(atom2) = rho_col(atom2) + rho_i_at_j
583 #else
584 rho(atom2) = rho(atom2) + rho_i_at_j
585 #endif
586 ! write(*,*) atom1,atom2,r,rho_i_at_j
587 endif
588
589 if (r.lt.EAMList%EAMParams(myid_atom2)%eam_rcut) then
590 call eam_splint(EAMList%EAMParams(myid_atom2)%eam_nr, &
591 EAMList%EAMParams(myid_atom2)%eam_rvals, &
592 EAMList%EAMParams(myid_atom2)%eam_rho_r, &
593 EAMList%EAMParams(myid_atom2)%eam_rho_r_pp, &
594 r, rho_j_at_i,drho,d2rho)
595
596
597
598
599 #ifdef IS_MPI
600 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
601 #else
602 rho(atom1) = rho(atom1) + rho_j_at_i
603 #endif
604 endif
605
606
607
608
609
610
611 end subroutine calc_eam_prepair_rho
612
613
614
615
616 !! Calculate the functional F(rho) for all local atoms
617 subroutine calc_eam_preforce_Frho(nlocal,pot)
618 integer :: nlocal
619 real(kind=dp) :: pot
620 integer :: i,j
621 integer :: atom
622 real(kind=dp) :: U,U1,U2
623 integer :: atype1
624 integer :: me
625 integer :: n_rho_points
626
627
628 cleanme = .true.
629 !! Scatter the electron density from pre-pair calculation back to local atoms
630 #ifdef IS_MPI
631 call scatter(rho_row,rho,plan_atom_row,eam_err)
632 if (eam_err /= 0 ) then
633 write(errMsg,*) " Error scattering rho_row into rho"
634 call handleError(RoutineName,errMesg)
635 endif
636 call scatter(rho_col,rho_tmp,plan_atom_col,eam_err)
637 if (eam_err /= 0 ) then
638 write(errMsg,*) " Error scattering rho_col into rho"
639 call handleError(RoutineName,errMesg)
640 endif
641
642 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
643 #endif
644
645
646
647 !! Calculate F(rho) and derivative
648 do atom = 1, nlocal
649 me = atid(atom)
650 n_rho_points = EAMList%EAMParams(me)%eam_nrho
651 ! Check to see that the density is not greater than the larges rho we have calculated
652 if (rho(atom) < EAMList%EAMParams(me)%eam_rhovals(n_rho_points)) then
653 call eam_splint(n_rho_points, &
654 EAMList%EAMParams(me)%eam_rhovals, &
655 EAMList%EAMParams(me)%eam_f_rho, &
656 EAMList%EAMParams(me)%eam_f_rho_pp, &
657 rho(atom), & ! Actual Rho
658 u, u1, u2)
659 else
660 ! Calculate F(rho with the largest available rho value
661 call eam_splint(n_rho_points, &
662 EAMList%EAMParams(me)%eam_rhovals, &
663 EAMList%EAMParams(me)%eam_f_rho, &
664 EAMList%EAMParams(me)%eam_f_rho_pp, &
665 EAMList%EAMParams(me)%eam_rhovals(n_rho_points), & ! Largest rho
666 u,u1,u2)
667 end if
668
669
670 frho(atom) = u
671 dfrhodrho(atom) = u1
672 d2frhodrhodrho(atom) = u2
673 pot = pot + u
674
675 enddo
676
677
678
679 #ifdef IS_MPI
680 !! communicate f(rho) and derivatives back into row and column arrays
681 call gather(frho,frho_row,plan_atom_row, eam_err)
682 if (eam_err /= 0) then
683 call handleError("cal_eam_forces()","MPI gather frho_row failure")
684 endif
685 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, eam_err)
686 if (eam_err /= 0) then
687 call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure")
688 endif
689 call gather(frho,frho_col,plan_atom_col, eam_err)
690 if (eam_err /= 0) then
691 call handleError("cal_eam_forces()","MPI gather frho_col failure")
692 endif
693 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, eam_err)
694 if (eam_err /= 0) then
695 call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure")
696 endif
697
698
699
700
701
702 if (nmflag) then
703 call gather(d2frhodrhodrho,d2frhodrhodrho_row,plan_atom_row)
704 call gather(d2frhodrhodrho,d2frhodrhodrho_col,plan_atom_col)
705 endif
706 #endif
707
708
709 end subroutine calc_eam_preforce_Frho
710
711
712
713
714 !! Does EAM pairwise Force calculation.
715 subroutine do_eam_pair(atom1, atom2, d, rij, r2, sw, vpair, fpair, &
716 pot, f, do_pot)
717 !Arguments
718 integer, intent(in) :: atom1, atom2
719 real( kind = dp ), intent(in) :: rij, r2
720 real( kind = dp ) :: pot, sw, vpair
721 real( kind = dp ), dimension(3,nLocal) :: f
722 real( kind = dp ), intent(in), dimension(3) :: d
723 real( kind = dp ), intent(inout), dimension(3) :: fpair
724
725 logical, intent(in) :: do_pot
726
727 real( kind = dp ) :: drdx,drdy,drdz
728 real( kind = dp ) :: d2
729 real( kind = dp ) :: phab,pha,dvpdr,d2vpdrdr
730 real( kind = dp ) :: rha,drha,d2rha, dpha
731 real( kind = dp ) :: rhb,drhb,d2rhb, dphb
732 real( kind = dp ) :: dudr
733 real( kind = dp ) :: rci,rcj
734 real( kind = dp ) :: drhoidr,drhojdr
735 real( kind = dp ) :: d2rhoidrdr
736 real( kind = dp ) :: d2rhojdrdr
737 real( kind = dp ) :: Fx,Fy,Fz
738 real( kind = dp ) :: r,d2pha,phb,d2phb
739
740 integer :: id1,id2
741 integer :: mytype_atom1
742 integer :: mytype_atom2
743
744 !Local Variables
745
746 ! write(*,*) "Frho: ", Frho(atom1)
747
748 phab = 0.0E0_DP
749 dvpdr = 0.0E0_DP
750 d2vpdrdr = 0.0E0_DP
751
752 if (rij .lt. EAM_rcut) then
753
754 #ifdef IS_MPI
755 mytype_atom1 = atid_row(atom1)
756 mytype_atom2 = atid_col(atom2)
757 #else
758 mytype_atom1 = atid(atom1)
759 mytype_atom2 = atid(atom2)
760 #endif
761 ! get cutoff for atom 1
762 rci = EAMList%EAMParams(mytype_atom1)%eam_rcut
763 ! get type specific cutoff for atom 2
764 rcj = EAMList%EAMParams(mytype_atom2)%eam_rcut
765
766 drdx = d(1)/rij
767 drdy = d(2)/rij
768 drdz = d(3)/rij
769
770 if (rij.lt.rci) then
771 call eam_splint(EAMList%EAMParams(mytype_atom1)%eam_nr, &
772 EAMList%EAMParams(mytype_atom1)%eam_rvals, &
773 EAMList%EAMParams(mytype_atom1)%eam_rho_r, &
774 EAMList%EAMParams(mytype_atom1)%eam_rho_r_pp, &
775 rij, rha,drha,d2rha)
776 !! Calculate Phi(r) for atom1.
777 call eam_splint(EAMList%EAMParams(mytype_atom1)%eam_nr, &
778 EAMList%EAMParams(mytype_atom1)%eam_rvals, &
779 EAMList%EAMParams(mytype_atom1)%eam_phi_r, &
780 EAMList%EAMParams(mytype_atom1)%eam_phi_r_pp, &
781 rij, pha,dpha,d2pha)
782 endif
783
784 if (rij.lt.rcj) then
785 ! Calculate rho,drho and d2rho for atom1
786 call eam_splint(EAMList%EAMParams(mytype_atom2)%eam_nr, &
787 EAMList%EAMParams(mytype_atom2)%eam_rvals, &
788 EAMList%EAMParams(mytype_atom2)%eam_rho_r, &
789 EAMList%EAMParams(mytype_atom2)%eam_rho_r_pp, &
790 rij, rhb,drhb,d2rhb)
791
792 !! Calculate Phi(r) for atom2.
793 call eam_splint(EAMList%EAMParams(mytype_atom2)%eam_nr, &
794 EAMList%EAMParams(mytype_atom2)%eam_rvals, &
795 EAMList%EAMParams(mytype_atom2)%eam_phi_r, &
796 EAMList%EAMParams(mytype_atom2)%eam_phi_r_pp, &
797 rij, phb,dphb,d2phb)
798 endif
799
800 if (rij.lt.rci) then
801 phab = phab + 0.5E0_DP*(rhb/rha)*pha
802 dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + &
803 pha*((drhb/rha) - (rhb*drha/rha/rha)))
804 d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rhb/rha)*d2pha + &
805 2.0E0_DP*dpha*((drhb/rha) - (rhb*drha/rha/rha)) + &
806 pha*((d2rhb/rha) - 2.0E0_DP*(drhb*drha/rha/rha) + &
807 (2.0E0_DP*rhb*drha*drha/rha/rha/rha) - (rhb*d2rha/rha/rha)))
808 endif
809
810 if (rij.lt.rcj) then
811 phab = phab + 0.5E0_DP*(rha/rhb)*phb
812 dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + &
813 phb*((drha/rhb) - (rha*drhb/rhb/rhb)))
814 d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rha/rhb)*d2phb + &
815 2.0E0_DP*dphb*((drha/rhb) - (rha*drhb/rhb/rhb)) + &
816 phb*((d2rha/rhb) - 2.0E0_DP*(drha*drhb/rhb/rhb) + &
817 (2.0E0_DP*rha*drhb*drhb/rhb/rhb/rhb) - (rha*d2rhb/rhb/rhb)))
818 endif
819
820 drhoidr = drha
821 drhojdr = drhb
822
823 d2rhoidrdr = d2rha
824 d2rhojdrdr = d2rhb
825
826
827 #ifdef IS_MPI
828 dudr = drhojdr*dfrhodrho_row(atom1)+drhoidr*dfrhodrho_col(atom2) &
829 + dvpdr
830
831 #else
832 dudr = drhojdr*dfrhodrho(atom1)+drhoidr*dfrhodrho(atom2) &
833 + dvpdr
834 ! write(*,*) "Atom1,Atom2, dfrhodrho(atom1) dfrhodrho(atom2): ", atom1,atom2,dfrhodrho(atom1),dfrhodrho(atom2)
835 #endif
836
837 fx = dudr * drdx
838 fy = dudr * drdy
839 fz = dudr * drdz
840
841
842 #ifdef IS_MPI
843 if (do_pot) then
844 pot_Row(atom1) = pot_Row(atom1) + phab*0.5
845 pot_Col(atom2) = pot_Col(atom2) + phab*0.5
846 end if
847
848 f_Row(1,atom1) = f_Row(1,atom1) + fx
849 f_Row(2,atom1) = f_Row(2,atom1) + fy
850 f_Row(3,atom1) = f_Row(3,atom1) + fz
851
852 f_Col(1,atom2) = f_Col(1,atom2) - fx
853 f_Col(2,atom2) = f_Col(2,atom2) - fy
854 f_Col(3,atom2) = f_Col(3,atom2) - fz
855 #else
856
857 if(do_pot) then
858 pot = pot + phab
859 end if
860
861 f(1,atom1) = f(1,atom1) + fx
862 f(2,atom1) = f(2,atom1) + fy
863 f(3,atom1) = f(3,atom1) + fz
864
865 f(1,atom2) = f(1,atom2) - fx
866 f(2,atom2) = f(2,atom2) - fy
867 f(3,atom2) = f(3,atom2) - fz
868 #endif
869
870 vpair = vpair + phab
871 #ifdef IS_MPI
872 id1 = AtomRowToGlobal(atom1)
873 id2 = AtomColToGlobal(atom2)
874 #else
875 id1 = atom1
876 id2 = atom2
877 #endif
878
879 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
880
881 fpair(1) = fpair(1) + fx
882 fpair(2) = fpair(2) + fy
883 fpair(3) = fpair(3) + fz
884
885 endif
886
887 if (nmflag) then
888
889 drhoidr = drha
890 drhojdr = drhb
891 d2rhoidrdr = d2rha
892 d2rhojdrdr = d2rhb
893
894 #ifdef IS_MPI
895 d2 = d2vpdrdr + &
896 d2rhoidrdr*dfrhodrho_col(atom2) + &
897 d2rhojdrdr*dfrhodrho_row(atom1) + &
898 drhoidr*drhoidr*d2frhodrhodrho_col(atom2) + &
899 drhojdr*drhojdr*d2frhodrhodrho_row(atom1)
900
901 #else
902
903 d2 = d2vpdrdr + &
904 d2rhoidrdr*dfrhodrho(atom2) + &
905 d2rhojdrdr*dfrhodrho(atom1) + &
906 drhoidr*drhoidr*d2frhodrhodrho(atom2) + &
907 drhojdr*drhojdr*d2frhodrhodrho(atom1)
908 #endif
909 end if
910
911 endif
912 end subroutine do_eam_pair
913
914
915 subroutine eam_splint(nx, xa, ya, yppa, x, y, dy, d2y)
916
917 integer :: atype, nx, j
918 real( kind = DP ), dimension(:) :: xa
919 real( kind = DP ), dimension(:) :: ya
920 real( kind = DP ), dimension(:) :: yppa
921 real( kind = DP ) :: x, y
922 real( kind = DP ) :: dy, d2y
923 real( kind = DP ) :: del, h, a, b, c, d
924 integer :: pp_arraySize
925
926
927 ! this spline code assumes that the x points are equally spaced
928 ! do not attempt to use this code if they are not.
929
930
931 ! find the closest point with a value below our own:
932 j = FLOOR(real((nx-1),kind=dp) * (x - xa(1)) / (xa(nx) - xa(1))) + 1
933
934 ! check to make sure we're inside the spline range:
935 if ((j.gt.nx).or.(j.lt.1)) then
936 write(errMSG,*) "EAM_splint: x is outside bounds of spline: ",x,j
937 call handleError(routineName,errMSG)
938 endif
939 ! check to make sure we haven't screwed up the calculation of j:
940 if ((x.lt.xa(j)).or.(x.gt.xa(j+1))) then
941 if (j.ne.nx) then
942 write(errMSG,*) "EAM_splint:",x," x is outside bounding range"
943 call handleError(routineName,errMSG)
944 endif
945 endif
946
947 del = xa(j+1) - x
948 h = xa(j+1) - xa(j)
949
950 a = del / h
951 b = 1.0E0_DP - a
952 c = a*(a*a - 1.0E0_DP)*h*h/6.0E0_DP
953 d = b*(b*b - 1.0E0_DP)*h*h/6.0E0_DP
954
955 y = a*ya(j) + b*ya(j+1) + c*yppa(j) + d*yppa(j+1)
956
957 dy = (ya(j+1)-ya(j))/h &
958 - (3.0E0_DP*a*a - 1.0E0_DP)*h*yppa(j)/6.0E0_DP &
959 + (3.0E0_DP*b*b - 1.0E0_DP)*h*yppa(j+1)/6.0E0_DP
960
961
962 d2y = a*yppa(j) + b*yppa(j+1)
963
964
965 end subroutine eam_splint
966
967
968 subroutine eam_spline(nx, xa, ya, yppa, yp1, ypn, boundary)
969
970
971 ! yp1 and ypn are the first derivatives of y at the two endpoints
972 ! if boundary is 'L' the lower derivative is used
973 ! if boundary is 'U' the upper derivative is used
974 ! if boundary is 'B' then both derivatives are used
975 ! if boundary is anything else, then both derivatives are assumed to be 0
976
977 integer :: nx, i, k, max_array_size
978
979 real( kind = DP ), dimension(:) :: xa
980 real( kind = DP ), dimension(:) :: ya
981 real( kind = DP ), dimension(:) :: yppa
982 real( kind = DP ), dimension(size(xa)) :: u
983 real( kind = DP ) :: yp1,ypn,un,qn,sig,p
984 character(len=*) :: boundary
985
986 ! make sure the sizes match
987 if ((nx /= size(xa)) .or. (nx /= size(ya))) then
988 call handleWarning("EAM_SPLINE","Array size mismatch")
989 end if
990
991 if ((boundary.eq.'l').or.(boundary.eq.'L').or. &
992 (boundary.eq.'b').or.(boundary.eq.'B')) then
993 yppa(1) = -0.5E0_DP
994 u(1) = (3.0E0_DP/(xa(2)-xa(1)))*((ya(2)-&
995 ya(1))/(xa(2)-xa(1))-yp1)
996 else
997 yppa(1) = 0.0E0_DP
998 u(1) = 0.0E0_DP
999 endif
1000
1001 do i = 2, nx - 1
1002 sig = (xa(i) - xa(i-1)) / (xa(i+1) - xa(i-1))
1003 p = sig * yppa(i-1) + 2.0E0_DP
1004 yppa(i) = (sig - 1.0E0_DP) / p
1005 u(i) = (6.0E0_DP*((ya(i+1)-ya(i))/(xa(i+1)-xa(i)) - &
1006 (ya(i)-ya(i-1))/(xa(i)-xa(i-1)))/ &
1007 (xa(i+1)-xa(i-1)) - sig * u(i-1))/p
1008 enddo
1009
1010 if ((boundary.eq.'u').or.(boundary.eq.'U').or. &
1011 (boundary.eq.'b').or.(boundary.eq.'B')) then
1012 qn = 0.5E0_DP
1013 un = (3.0E0_DP/(xa(nx)-xa(nx-1)))* &
1014 (ypn-(ya(nx)-ya(nx-1))/(xa(nx)-xa(nx-1)))
1015 else
1016 qn = 0.0E0_DP
1017 un = 0.0E0_DP
1018 endif
1019
1020 yppa(nx)=(un-qn*u(nx-1))/(qn*yppa(nx-1)+1.0E0_DP)
1021
1022 do k = nx-1, 1, -1
1023 yppa(k)=yppa(k)*yppa(k+1)+u(k)
1024 enddo
1025
1026 end subroutine eam_spline
1027
1028 end module eam