ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/suttonchen.F90
Revision: 2401
Committed: Tue Nov 1 18:29:57 2005 UTC (18 years, 8 months ago) by chuckv
File size: 32219 byte(s)
Log Message:
In process of adding sutton-chen forcefield.

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