ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/DarkSide/eam.F90
Revision: 2169
Committed: Tue Apr 12 17:12:27 2005 UTC (19 years, 2 months ago) by chuckv
File size: 32610 byte(s)
Log Message:
Updates to deallocate object in fortran.

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