ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/eam.F90
Revision: 2278
Committed: Fri Aug 26 22:39:25 2005 UTC (18 years, 10 months ago) by chrisfen
File size: 32664 byte(s)
Log Message:
updated getEAMCut

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