ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/eam.F90
Revision: 1930
Committed: Wed Jan 12 22:41:40 2005 UTC (19 years, 5 months ago) by gezelter
File size: 32713 byte(s)
Log Message:
merging new_design branch into OOPSE-2.0

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