ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/nano_mpi/src/force_eam_module.F90
Revision: 4
Committed: Mon Jun 10 17:18:36 2002 UTC (22 years ago) by chuckv
File size: 37942 byte(s)
Log Message:
Import Root

File Contents

# Content
1 module eam_module
2 use definitions, ONLY : DP,ndim,machdep_lnblnk
3 use parameter
4 use simulation
5 use second_deriv
6 use status, ONLY: error,info
7 use force_utilities, ONLY : wrap,check,save_nlist
8 #ifdef MPI
9 use mpi_module
10 use mpi_constants, ONLY: mpi_wtime
11 #endif
12
13
14
15
16 !! standard eam stuff
17 integer :: n_eam_atypes
18 integer, allocatable, dimension(:) :: eam_atype
19 real( kind = DP ), allocatable, dimension(:) :: eam_dr
20 integer, allocatable, dimension(:) :: eam_nr
21 integer, allocatable, dimension(:) :: eam_nrho
22 real( kind = DP ), allocatable, dimension(:) :: eam_lattice
23 real( kind = DP ), allocatable, dimension(:) :: eam_drho
24 integer , allocatable, dimension(:) :: eam_atype_map
25 real( kind = DP ), allocatable, dimension(:,:) :: eam_rvals
26 real( kind = DP ), allocatable, dimension(:,:) :: eam_rhovals
27 real( kind = DP ), allocatable, dimension(:) :: eam_rcut
28 real( kind = DP ), allocatable, dimension(:,:) :: eam_F_rho
29 real( kind = DP ), allocatable, dimension(:,:) :: eam_Z_r
30 real( kind = DP ), allocatable, dimension(:,:) :: eam_rho_r
31 real( kind = DP ), allocatable, dimension(:,:) :: eam_phi_r
32 real( kind = DP ), allocatable, dimension(:,:) :: eam_F_rho_pp
33 real( kind = DP ), allocatable, dimension(:,:) :: eam_Z_r_pp
34 real( kind = DP ), allocatable, dimension(:,:) :: eam_rho_r_pp
35 real( kind = DP ), allocatable, dimension(:,:) :: eam_phi_r_pp
36
37
38 real( kind = DP ), private :: time0,time1,time2,time3
39
40 integer, private :: eam_err
41
42 private :: mass_weight
43 private :: allocate_eam_atype,allocate_eam_module,deallocate_eam_module
44 private :: read_eam_pot, get_eam_sizes
45
46
47 contains
48
49
50 subroutine allocate_eam_atype(n_size_atype)
51 integer, intent(in) :: n_size_atype
52
53 allocate(eam_atype(n_size_atype))
54 allocate(eam_drho(n_size_atype))
55 allocate(eam_dr(n_size_atype))
56 allocate(eam_nr(n_size_atype))
57 allocate(eam_nrho(n_size_atype))
58 allocate(eam_lattice(n_size_atype))
59 allocate(eam_rcut(n_size_atype))
60
61 end subroutine allocate_eam_atype
62
63 subroutine allocate_eam_module(n_size_atype,n_eam_points)
64 integer, intent(in) :: n_eam_points
65 integer, intent(in) :: n_size_atype
66
67 allocate(eam_rvals(n_eam_points,n_size_atype))
68 allocate(eam_rhovals(n_eam_points,n_size_atype))
69 allocate(eam_F_rho(n_eam_points,n_size_atype))
70 allocate(eam_Z_r(n_eam_points,n_size_atype))
71 allocate(eam_rho_r(n_eam_points,n_size_atype))
72 allocate(eam_phi_r(n_eam_points,n_size_atype))
73 allocate(eam_F_rho_pp(n_eam_points,n_size_atype))
74 allocate(eam_Z_r_pp(n_eam_points,n_size_atype))
75 allocate(eam_rho_r_pp(n_eam_points,n_size_atype))
76 allocate(eam_phi_r_pp(n_eam_points,n_size_atype))
77
78 end subroutine allocate_eam_module
79
80 subroutine deallocate_eam_module()
81
82 deallocate(eam_atype)
83 deallocate(eam_drho)
84 deallocate(eam_dr)
85 deallocate(eam_nr)
86 deallocate(eam_nrho)
87 deallocate(eam_lattice)
88 deallocate(eam_atype_map)
89 deallocate(eam_rvals)
90 deallocate(eam_rhovals)
91 deallocate(eam_rcut)
92 deallocate(eam_Z_r)
93 deallocate(eam_rho_r)
94 deallocate(eam_phi_r)
95 deallocate(eam_F_rho_pp)
96 deallocate(eam_Z_r_pp)
97 deallocate(eam_rho_r_pp)
98 deallocate(eam_phi_r_pp)
99
100 end subroutine deallocate_eam_module
101
102 subroutine calc_eam_dens(update_nlist)
103
104 ! include 'headers/sizes.h'
105
106 real( kind = DP ) :: ptmp, rho_i_at_j,rho_j_at_i
107 real(kind=16) :: ptmp1, ptmp2
108 real(kind=16) :: this_rho, rho_total
109 integer :: i, j, atype1, atype2, nlist, jbeg, jend, jnab
110 integer :: tag_i,tag_j
111 real( kind = DP ) :: rxi, ryi, rzi, rxij, ryij, rzij, rijsq
112 real( kind = DP ) :: r, drho, d2rho
113 logical, intent(inout) :: update_nlist
114 integer :: j_start
115
116
117 #ifndef MPI
118 integer :: nrow
119 integer :: ncol
120
121 nrow = natoms - 1
122 ncol = natoms
123 #endif
124
125 this_rho = 0.0E0_DP
126 rho_total = 0.0E0_DP
127 #ifdef MPI
128 !Initialize timing info for MPI
129 time0 = 0.0E0_DP
130 time0 = mpi_wtime()
131
132
133 rho_row = 0.0E0_DP
134 rho_col = 0.0E0_DP
135 ! e_row = 0.0E0_DP
136 ! e_col = 0.0E0_DP
137
138 ! Initialize start of j index
139 j_start = 1
140
141 !distribute positions to row and col members
142 time1 = mpi_wtime()
143 call gather(q,q_row,plan_row3,eam_err)
144 if (eam_err /= 0) then
145 call error("calc_eam_dens()","MPI gather q_row failure")
146 end if
147
148
149
150 call gather(q,q_col,plan_col3,eam_err)
151 if (eam_err /= 0) then
152 call error("calc_eam_dens()","MPI gather q_col failure")
153 end if
154
155
156 time2 = mpi_wtime()
157 comm_time = comm_time + time2 - time1
158 #else
159 call cpu_time(time0)
160 rho = 0.0E0_DP
161 #endif
162
163 ! call mpi_barrier(mpi_comm_world,mpi_err)
164 if (update_nlist) then
165
166 ! save current configuration, contruct neighbor list,
167 ! and calculate forces
168
169 call save_nlist()
170
171
172 nlist = 0
173
174 do i = 1, nrow ! For normal nrow = natoms - 1
175 point(i) = nlist + 1
176
177 #ifdef MPI
178 tag_i = tag_row(i)
179 rxi = q_row(1,i)
180 ryi = q_row(2,i)
181 rzi = q_row(3,i)
182
183 #else
184 j_start = i + 1
185 rxi = q(1,i)
186 ryi = q(2,i)
187 rzi = q(3,i)
188 #endif
189
190 ! call mpi_barrier(mpi_comm_world,mpi_err)
191 inner: do j = j_start, ncol !For normal j_start is i + 1
192 !For MPI j_start is 1
193 !For normal ncol = 1
194 #ifdef MPI
195 ! call mpi_barrier(mpi_comm_world,mpi_err)
196 tag_j = tag_col(j)
197 if (newtons_thrd) then
198 if (tag_i <= tag_j) then
199 if (mod(tag_i + tag_j,2) == 0) cycle inner
200 else
201 if (mod(tag_i + tag_j,2) == 1) cycle inner
202 endif
203
204 endif
205
206
207
208 rxij = wrap(rxi - q_col(1,j), 1)
209 ryij = wrap(ryi - q_col(2,j), 2)
210 rzij = wrap(rzi - q_col(3,j), 3)
211
212 #else
213 rxij = wrap(rxi - q(1,j), 1)
214 ryij = wrap(ryi - q(2,j), 2)
215 rzij = wrap(rzi - q(3,j), 3)
216 #endif
217
218 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
219
220 !!$!_________________________remove me____________________________
221 !!$ r = dsqrt(rijsq)
222 !!$#ifdef MPI
223 !!$ if (tag_i == 1) then
224 !!$ if (node == 0) then
225 !!$ write(unit=*,fmt="(2i4,3es30.16)") tag_i,tag_j,rxij,ryij,rzij
226 !!$ endif
227 !!$ else if (tag_j == 1) then
228 !!$ if (node == 0) then
229 !!$ write(unit=*,fmt="(2i4,3es30.16)") tag_j,tag_i,&
230 !!$ rxij,ryij,rzij
231 !!$ endif
232 !!$ endif
233 !!$#else
234 !!$ if (i == 1) then
235 !!$ write(unit=*,fmt="(2i4,3es30.16)") i,j, &
236 !!$ rxij,ryij,rzij
237 !!$ endif
238 !!$
239 !!$
240 !!$#endif
241 !!$
242 !!$!----------------------end remove me
243 #ifdef MPI
244
245
246 if (rijsq <= rlstsq .AND. &
247 tag_j /= tag_i) then
248
249
250 #else
251 if (rijsq < rlstsq) then
252 #endif
253 nlist = nlist + 1
254 list(nlist) = j
255
256
257 if (rijsq <= rcutsq) then
258
259 r = dsqrt(rijsq)
260
261
262 ! find identities for i
263 #ifdef MPI
264 atype1 = ident_row(i)
265
266 #else
267 atype1 = ident(i)
268 #endif
269
270 ! find rho for atype1 - i
271 call calc_eam_rho(r, rho_i_at_j, drho, d2rho, atype1)
272
273
274 ! density at site j depends on type of atom at site i
275 #ifdef MPI
276
277 rho_col(j) = rho_col(j) + rho_i_at_j
278 ptmp1 = rho_i_at_j
279 #else
280 rho(j) = rho(j) + rho_i_at_j
281 ptmp1 = rho_i_at_j
282 #endif
283
284 ! find identities for j
285 #ifdef MPI
286 atype2 = ident_col(j)
287 #else
288 atype2 = ident(j)
289 #endif
290
291
292 ! find rho for atype2 - j
293 call calc_eam_rho(r, rho_j_at_i, drho, d2rho, atype2)
294
295
296 ! density at site i depends on type of atom at site j
297 #ifdef MPI
298
299 rho_row(i) = rho_row(i) + rho_j_at_i
300 ptmp2 = rho_j_at_i
301 #else
302 rho(i) = rho(i) + rho_j_at_i
303 ptmp2 = rho_j_at_i
304 #endif
305
306 !!$#ifdef MPI
307 !!$ if (tag_i == 1) then
308 !!$ if (node == 0) then
309 !!$ write(unit=*,fmt="(2i4,4es30.16)") tag_i,tag_j,r,rijsq,&
310 !!$ rho_i_at_j,rho_j_at_i
311 !!$ this_rho = this_rho + rho_j_at_i
312 !!$ endif
313 !!$ else if (tag_j == 1) then
314 !!$ if (node == 0) then
315 !!$ write(unit=*,fmt="(2i4,4es30.16)") tag_j,tag_i,r,rijsq,&
316 !!$ rho_i_at_j,rho_j_at_i
317 !!$ this_rho = this_rho + rho_i_at_j
318 !!$ endif
319 !!$ endif
320 !!$#else
321 !!$ if (i == 1) then
322 !!$ write(unit=*,fmt="(2i4,4es30.16)") i,j,r,rijsq, &
323 !!$ rho_i_at_j, rho_j_at_i
324 !!$ this_rho = this_rho + rho_j_at_i
325 !!$ endif
326 !!$
327 !!$
328 !!$#endif
329
330 endif
331 endif
332
333 enddo inner
334 end do
335
336
337
338 #ifdef MPI
339 point(nrow+1) = nlist + 1
340 #else
341 point(natoms) = nlist + 1
342 #endif
343
344 else
345 ! use the list to find the neighbors
346
347 do i = 1, nrow !nrow for non-MPI is just natoms
348 JBEG = POINT(I)
349 JEND = POINT(I+1) - 1
350
351 ! check thiat molecule i has neighbors
352 if (jbeg <= jend) then
353 #ifdef MPI
354 rxi = q_row(1,i)
355 ryi = q_row(2,i)
356 rzi = q_row(3,i)
357 #else
358 rxi = q(1,i)
359 ryi = q(2,i)
360 rzi = q(3,i)
361 #endif
362
363 do jnab = jbeg, jend
364 j = list(jnab)
365 #ifdef MPI
366 rxij = wrap(rxi - q_col(1,j), 1)
367 ryij = wrap(ryi - q_col(2,j), 2)
368 rzij = wrap(rzi - q_col(3,j), 3)
369 #else
370 rxij = wrap(rxi - q(1,j), 1)
371 ryij = wrap(ryi - q(2,j), 2)
372 rzij = wrap(rzi - q(3,j), 3)
373 #endif
374
375 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
376
377 if (rijsq < rcutsq) then
378
379 r = dsqrt(rijsq)
380
381 ! Find identities for i and j
382 #ifdef MPI
383 atype1 = ident_row(i)
384 atype2 = ident_col(j)
385 #else
386 atype1 = ident(i)
387 atype2 = ident(j)
388 #endif
389
390 call calc_eam_rho(r, ptmp, drho, d2rho, atype1)
391
392 ! density at site j depends on type of atom at site i
393 #ifdef MPI
394 rho_col(j) = rho_col(j) + ptmp
395 #else
396 rho(j) = rho(j) + ptmp
397 #endif
398
399 call calc_eam_rho(r, ptmp, drho, d2rho, atype2)
400
401 ! density at site i depends on type of atom at site j
402 #ifdef MPI
403 rho_row(i) = rho_row(i) + ptmp
404 #else
405 rho(i) = rho(i) + ptmp
406 #endif
407 endif
408 enddo
409 endif
410 enddo
411 endif
412
413 #ifdef MPI
414
415 !! communicate densities
416 time1 = mpi_wtime()
417 call mpi_barrier(mpi_comm_world,mpi_err)
418 ! write(*,*) "This rho", this_rho
419 ! call mpi_allreduce(this_rho, rho_total,1,mpi_double_precision, &
420 ! mpi_sum,mpi_comm_world,mpi_err)
421 ! if (node == 0) then
422 ! write(*,'(a27,es30.16)') "Rho total for particle 1", rho_total
423 ! endif
424 call mpi_barrier(mpi_comm_world,mpi_err)
425 call scatter(rho_row,rho,plan_row,eam_err)
426 if (eam_err /= 0) then
427 call error("calc_eam_dens()","MPI scatter rho_row failure")
428 endif
429 ! if (node == 0 ) then
430 ! write(*,*) "Rho before col comm: ", rho(1)
431 ! endif
432
433 if (newtons_thrd) then
434 call scatter(rho_col,rho_tmp,plan_col,eam_err)
435 if (eam_err /= 0) then
436 call error("calc_eam_dens()","MPI scatter rho_col failure")
437 endif
438 ! if (node == 0 ) then
439 ! write(*,*) "Rho tmp col: ", rho_tmp(1)
440 ! endif
441 do i = 1, nlocal
442 rho(i) = rho(i) + rho_tmp(i)
443 end do
444 ! if (node == 0 ) then
445 ! write(*,*) "Rho after col comm: ", rho(1)
446 ! endif
447 endif
448 time2 = mpi_wtime()
449 comm_time = comm_time + time2 - time1
450 #else
451 ! write(*,*) this_rho
452 #endif
453
454 return
455 end subroutine calc_eam_dens
456
457 subroutine calc_eam_forces(nmflag,pot)
458
459
460 ! include 'headers/sizes.h'
461 ! include 'headers/fileio.h'
462
463 #ifdef MPI
464 ! real( kind = DP ), dimension(nlocal) :: frho
465 real( kind = DP ), dimension(nlocal) :: dfrhodrho
466 real( kind = DP ), dimension(nlocal) :: d2frhodrhodrho
467 real( kind = DP ), dimension(3,ncol) :: efr
468 #else
469 real( kind = DP ), dimension(natoms) :: frho
470 real( kind = DP ), dimension(natoms) :: dfrhodrho
471 real( kind = DP ), dimension(natoms) :: d2frhodrhodrho
472 real( kind = DP ), dimension(3,natoms) :: efr
473 #endif
474
475
476 real( kind = DP ), intent(out), optional :: pot
477 real( kind = DP ) :: vptmp, dudr, ftmp
478 real( kind = DP ) :: u, u1, u2, phab, rci, rcj
479 real( kind = DP ) :: rha, drha, d2rha, pha, dpha, d2pha
480 real( kind = DP ) :: rhb, drhb, d2rhb, phb, dphb, d2phb
481 real( kind = DP ) :: drhoidr, drhojdr, d2rhoidrdr, d2rhojdrdr
482 real( kind = DP ) :: dvpdr, drdx1, d2vpdrdr, d2
483 real( kind = DP ) :: kt1, kt2, kt3, ktmp
484
485 real( kind = DP ) :: col_pot_total
486
487 integer :: i, j, dim, atype1, atype2, idim, jdim, dim2, idim2, jdim2
488 integer :: jbeg, jend, jnab
489 real( kind = DP ) :: rxij, ryij, rzij, rxi, ryi, rzi, rijsq, r
490
491 integer :: nlist
492 logical, intent(in) :: nmflag
493 logical :: do_pot
494
495 ! MPI variables and arrays.
496 #ifdef MPI
497 real( kind = DP ), dimension(nrow) :: frho_row
498 real( kind = DP ), dimension(ncol) :: frho_col
499 real( kind = DP ), dimension(nrow) :: dfrhodrho_row
500 real( kind = DP ), dimension(nrow) :: d2frhodrhodrho_row
501 real( kind = DP ), dimension(ncol) :: dfrhodrho_col
502 real( kind = DP ), dimension(ncol) :: d2frhodrhodrho_col
503
504 real( kind = DP ) :: pot_local, pot_phi_row, pot_Frho, pot_phi, pot_row
505 ! real( kind = DP ) :: pot1, pot2
506 ! normal variables
507 #else
508 integer :: nrow
509 integer :: ncol
510
511 nrow = natoms - 1
512 ncol = natoms
513 #endif
514
515
516
517 ! figure out if we should do pot.
518 do_pot = .false.
519 if (present(pot)) do_pot = .true.
520
521 if (do_pot) pot = 0.0E0_DP
522
523 #ifndef MPI
524 f = 0.0E0_DP
525 e = 0.0E0_DP
526 #else
527
528 f_row = 0.0E0_DP
529 f_col = 0.0E0_DP
530
531 pot_phi_row = 0.0E0_DP
532 pot_phi = 0.0E0_DP
533 pot_Frho = 0.0E0_DP
534 pot_local = 0.0E0_DP
535 pot_row = 0.0E0_DP
536
537 e_row = 0.0E0_DP
538 e_col = 0.0E0_DP
539 e_tmp = 0.0E0_DP
540 #endif
541
542 do i = 1, nlocal
543 atype1 = ident(i)
544
545 call calc_eam_frho(rho(i), u, u1, u2, atype1)
546 frho(i) = u
547 dfrhodrho(i) = u1
548 d2frhodrhodrho(i) = u2
549 ! pot_local = pot_local + u
550 #ifndef MPI
551 if (do_pot) pot = pot + u
552 #endif
553 enddo
554
555 #ifdef MPI
556 time1 = mpi_wtime()
557 !! communicate f(rho) and derivatives
558 call gather(frho,frho_row,plan_row, eam_err)
559 if (eam_err /= 0) then
560 call error("cal_eam_forces()","MPI gather frho_row failure")
561 endif
562 call gather(dfrhodrho,dfrhodrho_row,plan_row, eam_err)
563 if (eam_err /= 0) then
564 call error("cal_eam_forces()","MPI gather dfrhodrho_row failure")
565 endif
566 call gather(frho,frho_col,plan_col, eam_err)
567 if (eam_err /= 0) then
568 call error("cal_eam_forces()","MPI gather frho_col failure")
569 endif
570 call gather(dfrhodrho,dfrhodrho_col,plan_col, eam_err)
571 if (eam_err /= 0) then
572 call error("cal_eam_forces()","MPI gather dfrhodrho_col failure")
573 endif
574
575 if (nmflag) then
576 call gather(d2frhodrhodrho,d2frhodrhodrho_row,plan_row)
577 call gather(d2frhodrhodrho,d2frhodrhodrho_col,plan_col)
578 endif
579 time2 = mpi_wtime()
580 comm_time = comm_time + time2 - time1
581 #endif
582
583 do i = 1, nrow ! for normal nrow = natoms - 1
584 JBEG = POINT(i)
585 JEND = POINT(i+1) - 1
586
587 ! check thiat molecule i has neighbors
588 if (jbeg .le. jend) then
589 #ifdef MPI
590 atype1 = ident_row(i)
591 rxi = q_row(1,i)
592 ryi = q_row(2,i)
593 rzi = q_row(3,i)
594 #else
595 atype1 = ident(i)
596 rxi = q(1,i)
597 ryi = q(2,i)
598 rzi = q(3,i)
599 #endif
600
601 do jnab = jbeg, jend
602 j = list(jnab)
603 #ifdef MPI
604 rxij = wrap(rxi - q_col(1,j), 1)
605 ryij = wrap(ryi - q_col(2,j), 2)
606 rzij = wrap(rzi - q_col(3,j), 3)
607 #else
608 rxij = wrap(rxi - q(1,j), 1)
609 ryij = wrap(ryi - q(2,j), 2)
610 rzij = wrap(rzi - q(3,j), 3)
611 #endif
612
613 rijsq = rxij*rxij + ryij*ryij + rzij*rzij
614
615 if (rijsq .lt. rcutsq) then
616
617 r = dsqrt(rijsq)
618 efr(1,j) = -rxij
619 efr(2,j) = -ryij
620 efr(3,j) = -rzij
621
622
623 call calc_eam_rho(r, rha, drha, d2rha, atype1)
624 call calc_eam_phi(r, pha, dpha, d2pha, atype1)
625 rci = eam_rcut(eam_atype_map(atype1))
626 #ifdef MPI
627 atype2 = ident_col(j)
628 #else
629 atype2 = ident(j)
630 #endif
631
632 call calc_eam_rho(r, rhb, drhb, d2rhb, atype2)
633 call calc_eam_phi(r, phb, dphb, d2phb, atype2)
634 rcj = eam_rcut(eam_atype_map(atype2))
635
636 phab = 0.0E0_DP
637 dvpdr = 0.0E0_DP
638 d2vpdrdr = 0.0E0_DP
639
640 if (r.lt.rci) then
641 phab = phab + 0.5E0_DP*(rhb/rha)*pha
642 dvpdr = dvpdr + 0.5E0_DP*((rhb/rha)*dpha + &
643 pha*((drhb/rha) - (rhb*drha/rha/rha)))
644 d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rhb/rha)*d2pha + &
645 2.0E0_DP*dpha*((drhb/rha) - (rhb*drha/rha/rha)) + &
646 pha*((d2rhb/rha) - 2.0E0_DP*(drhb*drha/rha/rha) + &
647 (2.0E0_DP*rhb*drha*drha/rha/rha/rha) - (rhb*d2rha/rha/rha)))
648 endif
649
650
651 if (r.lt.rcj) then
652 phab = phab + 0.5E0_DP*(rha/rhb)*phb
653 dvpdr = dvpdr + 0.5E0_DP*((rha/rhb)*dphb + &
654 phb*((drha/rhb) - (rha*drhb/rhb/rhb)))
655 d2vpdrdr = d2vpdrdr + 0.5E0_DP*((rha/rhb)*d2phb + &
656 2.0E0_DP*dphb*((drha/rhb) - (rha*drhb/rhb/rhb)) + &
657 phb*((d2rha/rhb) - 2.0E0_DP*(drha*drhb/rhb/rhb) + &
658 (2.0E0_DP*rha*drhb*drhb/rhb/rhb/rhb) - (rha*d2rhb/rhb/rhb)))
659 endif
660
661
662 #ifdef MPI
663
664 e_row(i) = e_row(i) + phab*0.5
665 e_col(i) = e_col(i) + phab*0.5
666 #else
667 if (do_pot) pot = pot + phab
668 #endif
669
670 drhoidr = drha
671 drhojdr = drhb
672
673 d2rhoidrdr = d2rha
674 d2rhojdrdr = d2rhb
675 #ifdef MPI
676 dudr = drhojdr*dfrhodrho_row(i)+drhoidr*dfrhodrho_col(j) &
677 + dvpdr
678
679 if (nmflag) then
680 d2 = d2vpdrdr + &
681 d2rhoidrdr*dfrhodrho_col(j) + &
682 d2rhojdrdr*dfrhodrho_row(i) + &
683 drhoidr*drhoidr*d2frhodrhodrho_col(j) + &
684 drhojdr*drhojdr*d2frhodrhodrho_row(i)
685 endif
686 #else
687 dudr = drhojdr*dfrhodrho(i)+drhoidr*dfrhodrho(j) &
688 + dvpdr
689
690 d2 = d2vpdrdr + &
691 d2rhoidrdr*dfrhodrho(j) + &
692 d2rhojdrdr*dfrhodrho(i) + &
693 drhoidr*drhoidr*d2frhodrhodrho(j) + &
694 drhojdr*drhojdr*d2frhodrhodrho(i)
695 #endif
696
697
698 do dim = 1, 3
699
700 drdx1 = efr(dim,j) / r
701 ftmp = dudr * drdx1
702
703 #ifdef MPI
704 f_col(dim,j) = f_col(dim,j) - ftmp
705 f_row(dim,i) = f_row(dim,i) + ftmp
706 #else
707 f(dim,j) = f(dim,j) - ftmp
708 f(dim,i) = f(dim,i) + ftmp
709 #endif
710
711 if (nmflag) then
712 idim = 3 * (i-1) + dim
713 jdim = 3 * (j-1) + dim
714
715 do dim2 = 1, 3
716
717 kt1 = d2 * efr(dim,j) * efr(dim2,j)/r/r
718 kt2 = - dudr * efr(dim,j) * efr(dim2,j)/r/r/r
719
720 if (dim.eq.dim2) then
721 kt3 = dudr / r
722 else
723 kt3 = 0.0E0_DP
724 endif
725
726 ! The factor of 2 below is to compensate for
727 ! overcounting.
728 ! Mass weighting is done separately...
729
730 ktmp = (kt1+kt2+kt3)/2.0E0_DP
731 idim2 = 3 * (i-1) + dim2
732 jdim2 = 3 * (j-1) + dim2
733
734 d(idim, idim2) = d(idim,idim2) + ktmp
735 d(idim2, idim) = d(idim2,idim) + ktmp
736
737 d(idim, jdim2) = d(idim,jdim2) - ktmp
738 d(idim2, jdim) = d(idim2,jdim) - ktmp
739
740 d(jdim, idim2) = d(jdim,idim2) - ktmp
741 d(jdim2, idim) = d(jdim2,idim) - ktmp
742
743 d(jdim, jdim2) = d(jdim,jdim2) + ktmp
744 d(jdim2, jdim) = d(jdim2,jdim) + ktmp
745
746 enddo
747 endif
748 enddo
749
750 endif
751 enddo
752 endif
753
754 enddo
755
756
757 #ifdef MPI
758 time1 = mpi_wtime()
759 !!distribute forces
760 call scatter(f_row,f,plan_row3,eam_err)
761 if (eam_err /= 0) then
762 call error("calc_eam_forces()","MPI scatter f_row failure")
763 endif
764 if (newtons_thrd) then
765 call scatter(f_col,f_tmp,plan_col3,eam_err)
766 if (eam_err /= 0) then
767 call error("calc_eam_forces()","MPI scatter f_col failure")
768 endif
769 do i = 1,nlocal
770 do dim = 1,3
771 f(dim,i) = f(dim,i) + f_tmp(dim,i)
772 end do
773 end do
774 endif
775
776
777 if (do_pot) then
778 ! scatter/gather pot_row into the members of my column
779 call scatter(e_row,e_tmp,plan_row,eam_err)
780 if (eam_err /= 0) then
781 call error("calc_eam_forces()","MPI scatter e_row failure")
782 endif
783
784 ! scatter/gather pot_local into all other procs
785 ! add resultant to get total pot
786 do i = 1, nlocal
787 pot_local = pot_local + frho(i) + e_tmp(i)
788 enddo
789
790
791 if (newtons_thrd) then
792 e_tmp = 0.0E0_DP
793 call scatter(e_col,e_tmp,plan_col,eam_err)
794 if (eam_err /= 0) then
795 call error("calc_eam_forces()","MPI scatter e_col failure")
796 endif
797
798 do i = 1, nlocal
799 pot_local = pot_local + e_tmp(i)
800 enddo
801 endif
802
803 call mpi_reduce(pot_local,pot,1,mpi_double_precision, &
804 mpi_sum,0,mpi_comm_world,mpi_err)
805 if (mpi_err /= 0) then
806 call error("EAM_MODULE","MPI reduce pot_local error")
807 endif
808
809 endif
810
811 time2 = mpi_wtime()
812 comm_time = comm_time + time2 - time1
813 force_time = force_time + time2 - time0
814 #else
815 call cpu_time(time2)
816 force_time = force_time + time2 - time0
817 #endif
818
819
820 if (nmflag) then
821 call mass_weight()
822 endif
823
824 return
825 end subroutine calc_eam_forces
826
827 subroutine initialize_eam()
828 use model_module
829 use file_units, ONLY : next_unit
830
831
832
833 character(len=80) :: eam_pot_file
834 integer :: i, j, max_size, prev_max_size
835 integer :: number_rho, number_r
836 integer :: eam_unit
837 integer :: this_error
838 character(len=300) :: msg
839 integer, external :: nfiles
840 !for mpi
841
842
843 #ifdef MPI
844 if (node == 0) &
845 n_eam_atypes = nfiles(trim(eam_pot_dir)//char(0))
846
847 call mpi_bcast(n_eam_atypes,1,mpi_integer,0,mpi_comm_world,mpi_err)
848 if (n_eam_atypes == -1) then
849 call error("INITIALIZE_EAM","NO EAM potentials found!")
850 endif
851 write(msg,'(a5,i4,a12,i5,a14)') 'Node: ',node,' reading ...', &
852 n_eam_atypes, ' eam atom types'
853 call info('INITIALIZE_EAM', trim(msg))
854 #else
855 n_eam_atypes = nfiles(trim(eam_pot_dir)//char(0))
856 if (n_eam_atypes == -1) then
857 call error("INITIALIZE_EAM","NO EAM potentials found!")
858 endif
859
860 write(msg,'(a12,i5,a14)') ' Reading ...', &
861 n_eam_atypes, ' eam atom types'
862 call info('INITIALIZE_EAM', trim(msg))
863 #endif
864
865
866 call allocate_eam_atype(n_eam_atypes)
867
868
869
870 !! get largest number of data points for any potential
871 #ifdef MPI
872 if (node == 0) then
873 #endif
874 prev_max_size = 0
875 do i = 1, n_eam_atypes
876 call getfilename(i, eam_pot_file)
877 max_size = max(get_eam_sizes( &
878 trim(eam_pot_dir) // '/' // eam_pot_file), &
879 prev_max_size)
880 prev_max_size = max_size
881 end do
882 #ifdef MPI
883 end if
884
885
886 call mpi_bcast(max_size,1,mpi_integer,0,mpi_comm_world,mpi_err)
887 #endif
888
889 call allocate_eam_module(n_eam_atypes,max_size)
890 allocate(eam_atype_map(get_max_atype()))
891
892 #ifdef MPI
893 if (node == 0) then
894 #endif
895 do i = 1, n_eam_atypes
896 call getfilename(i, eam_pot_file)
897 call read_eam_pot(i,trim(eam_pot_dir) // '/' // eam_pot_file, &
898 this_error)
899
900 do j = 1, eam_nr(i)
901 eam_rvals(j,i) = dble(j-1)*eam_dr(i)
902 enddo
903
904 do j = 1, eam_nrho(i)
905 eam_rhovals(j,i) = dble(j-1)*eam_drho(i)
906 enddo
907
908 ! convert from eV to kcal / mol:
909 do j = 1, eam_nrho(i)
910 eam_F_rho(j,i) = eam_F_rho(j,i)*23.06054E0_DP
911 enddo
912
913 ! precompute the pair potential and get it into kcal / mol:
914 eam_phi_r(1,i) = 0.0E0_DP
915 do j = 2, eam_nr(i)
916 eam_phi_r(j,i) = (eam_Z_r(j,i)**2)/eam_rvals(j,i)
917 eam_phi_r(j,i) = eam_phi_r(j,i)*331.999296E0_DP
918 enddo
919
920 end do
921 #ifdef MPI
922 call info('INITIALIZE_EAM','NODE 0: Distributing spline arrays')
923 endif
924
925 call mpi_bcast(this_error,n_eam_atypes,mpi_integer,0, &
926 mpi_comm_world,mpi_err)
927 if (this_error /= 0) then
928 call error('INITIALIZE_EAM',"Cannot read eam files")
929 endif
930
931 call mpi_bcast(eam_atype,n_eam_atypes,mpi_integer,0, &
932 mpi_comm_world,mpi_err)
933
934 !! distribute values to cluster......
935 call mpi_bcast(eam_nr,n_eam_atypes,mpi_integer,&
936 0,mpi_comm_world,mpi_err)
937 call mpi_bcast(eam_nrho,n_eam_atypes,mpi_integer,&
938 0,mpi_comm_world,mpi_err)
939 call mpi_bcast(eam_rvals,n_eam_atypes*max_size,mpi_double_precision, &
940 0,mpi_comm_world,mpi_err)
941 call mpi_bcast(eam_rcut,n_eam_atypes,mpi_double_precision, &
942 0,mpi_comm_world,mpi_err)
943 call mpi_bcast(eam_rhovals,n_eam_atypes*max_size,mpi_double_precision, &
944 0,mpi_comm_world,mpi_err)
945
946 !! distribute arrays
947 call mpi_bcast(eam_rho_r,n_eam_atypes*max_size,mpi_double_precision, &
948 0,mpi_comm_world,mpi_err)
949 call mpi_bcast(eam_Z_r,n_eam_atypes*max_size,mpi_double_precision, &
950 0,mpi_comm_world,mpi_err)
951 call mpi_bcast(eam_F_rho,n_eam_atypes*max_size,mpi_double_precision, &
952 0,mpi_comm_world,mpi_err)
953 call mpi_bcast(eam_phi_r,n_eam_atypes*max_size,mpi_double_precision, &
954 0,mpi_comm_world,mpi_err)
955
956 #endif
957 call info('INITIALIZE_EAM', 'creating splines')
958
959 do i = 1, n_eam_atypes
960 number_r = eam_nr(i)
961 number_rho = eam_nrho(i)
962
963 call eam_spline(i, number_r, eam_rvals, eam_rho_r, eam_rho_r_pp, &
964 0.0E0_DP, 0.0E0_DP, 'N')
965 call eam_spline(i, number_r, eam_rvals, eam_Z_r, eam_Z_r_pp, &
966 0.0E0_DP, 0.0E0_DP, 'N')
967 call eam_spline(i, number_rho, eam_rhovals, eam_F_rho, eam_F_rho_pp, &
968 0.0E0_DP, 0.0E0_DP, 'N')
969 call eam_spline(i, number_r, eam_rvals, eam_phi_r, eam_phi_r_pp, &
970 0.0E0_DP, 0.0E0_DP, 'N')
971 enddo
972
973 do i = 1, n_eam_atypes
974 eam_atype_map(eam_atype(i)) = i
975 end do
976
977
978
979 call info('INITIALIZE_EAM','Done creating splines')
980
981 return
982 end subroutine initialize_eam
983
984 subroutine read_eam_pot(atype_index, eam_pot_file,error)
985 use model_module
986 use file_units, ONLY : next_unit
987
988
989 !! variables to store number of potential points
990
991 integer :: j
992 integer :: ierr
993 integer :: atype_index
994 integer :: junk
995 integer :: pot_unit
996 integer, intent(out), optional :: error
997 !! character format parameters for Dynamo files
998 character(len=*), parameter :: line_2 = "(i5,2d15.15)"
999 character(len=*), parameter :: line = "(i5,d24.16,i5,d24.16,d24.16)"
1000 ! character(len=*), parameter :: line_8 = "*"
1001 character(len=*), parameter :: potential_lines = "(5d24.16)"
1002 character(len=*), intent(in) :: eam_pot_file
1003 character(len=80) :: msg
1004 character(len=80) :: junkline
1005 real( kind = DP ) :: junk1, junk2
1006
1007
1008 error = 0
1009
1010 ! open potential file first
1011 pot_unit = next_unit()
1012 open (unit=pot_unit,file=eam_pot_file, &
1013 status="old",iostat=ierr,action="read")
1014 !! handle error if file cannot be opened....
1015 if (ierr > 0) then
1016 write(msg,*) "Error opening potential model file", trim(eam_pot_file)
1017 call info('read_eam_pot', msg)
1018 error = -1
1019 end if
1020
1021 !------------------> read top of file
1022 read(pot_unit,'(a80)') junkline ! first line is a comment line
1023 ! read line 2 atomic number, atomic mass and lattice constant
1024 read(pot_unit,line_2) eam_atype(atype_index), junk2, &
1025 eam_lattice(atype_index)
1026
1027
1028 write(msg,*) 'Found potential for: ', atype_identifier(eam_atype(atype_index))
1029 call info('READ_EAM_POT', msg)
1030
1031 ! read line 3
1032 read(pot_unit,*) &
1033 eam_nrho(atype_index), &
1034 eam_drho(atype_index), &
1035 eam_nr(atype_index), &
1036 eam_dr(atype_index), &
1037 eam_rcut(atype_index)
1038
1039 !-------------------> read potential points
1040
1041 ! read F(rho)
1042 read(pot_unit,potential_lines) &
1043 (eam_F_rho(j,atype_index),j=1,eam_nrho(atype_index))
1044
1045 !! read in Z(r)
1046 read(pot_unit,potential_lines) &
1047 (eam_Z_r(j,atype_index),j=1,eam_nr(atype_index))
1048
1049 !! read in rho(r)
1050 read(pot_unit,potential_lines) &
1051 (eam_rho_r(j,atype_index),j=1,eam_nr(atype_index))
1052
1053 close (pot_unit)
1054
1055 end subroutine read_eam_pot
1056
1057 function get_eam_sizes(eam_pot_file) result(max_size)
1058 use file_units, ONLY : next_unit
1059
1060 character(len=*) :: eam_pot_file
1061 character(len=*), parameter :: line_3 = "(i5,f24.16,i5,2f24.16)"
1062 integer :: j
1063 integer :: ierr
1064 integer :: atype_index
1065 integer :: machdep_lnblnk
1066 integer :: number_of_rho
1067 integer :: number_of_r
1068 integer :: max_size
1069 integer :: pot_unit
1070 real( kind = DP ) :: junk1,junk2,junk3
1071 character(len=80) :: msg,junk
1072 ! open potential file first
1073 pot_unit = next_unit()
1074 open (unit=pot_unit,file=eam_pot_file, &
1075 status="old",iostat=ierr,action="read")
1076 !! handle error if file cannot be opened....
1077 if (ierr > 0) then
1078 write(msg,*) "Error opening potential model file", trim(eam_pot_file)
1079 call error('read_eam_pot', msg)
1080 end if
1081
1082 read(pot_unit,*) ! first line is a comment line
1083 ! read line 2 atomic number, atomic mass and lattice constant
1084 read(pot_unit,*)
1085
1086 read(pot_unit,line_3) &
1087 number_of_rho, &
1088 junk1, &
1089 number_of_r, &
1090 junk2, &
1091 junk3
1092
1093 close(pot_unit)
1094
1095 max_size = max(number_of_rho,number_of_r)
1096
1097
1098 end function get_eam_sizes
1099
1100
1101
1102 subroutine eam_splint(atype, nx, xa, ya, yppa, x, y, dy, d2y)
1103
1104 ! include 'headers/sizes.h'
1105
1106 real( kind = DP ), dimension(:,:) :: xa
1107 real( kind = DP ), dimension(:,:) :: ya
1108 real( kind = DP ), dimension(:,:) :: yppa
1109 real( kind = DP ) :: x, y, dy, d2y
1110 real( kind = DP ) :: del, h, a, b, c, d
1111
1112
1113 integer atype, nx, j
1114
1115
1116 ! this spline code assumes that the x points are equally spaced
1117 ! do not attempt to use this code if they are not.
1118
1119
1120 ! find the closest point with a value below our own:
1121 j = FLOOR(dble(nx-1) * (x - xa(1,atype)) / (xa(nx,atype) - xa(1,atype))) + 1
1122
1123 ! check to make sure we're inside the spline range:
1124 if ((j.gt.nx).or.(j.lt.1)) call error('eam_splint', &
1125 'x is outside bounds of spline')
1126
1127 ! check to make sure we haven't screwed up the calculation of j:
1128 if ((x.lt.xa(j,atype)).or.(x.gt.xa(j+1,atype))) then
1129 if (j.ne.nx) then
1130 call error('eam_splint', &
1131 'x is outside bounding range')
1132 endif
1133 endif
1134
1135 del = xa(j+1,atype) - x
1136 h = xa(j+1,atype) - xa(j,atype)
1137
1138 a = del / h
1139 b = 1.0E0_DP - a
1140 c = a*(a*a - 1.0E0_DP)*h*h/6.0E0_DP
1141 d = b*(b*b - 1.0E0_DP)*h*h/6.0E0_DP
1142
1143 y = a*ya(j,atype) + b*ya(j+1,atype) + c*yppa(j,atype) + d*yppa(j+1,atype)
1144
1145 dy = (ya(j+1,atype)-ya(j,atype))/h &
1146 - (3.0E0_DP*a*a - 1.0E0_DP)*h*yppa(j,atype)/6.0E0_DP &
1147 + (3.0E0_DP*b*b - 1.0E0_DP)*h*yppa(j+1,atype)/6.0E0_DP
1148
1149 d2y = a*yppa(j,atype) + b*yppa(j+1,atype)
1150
1151 return
1152 end subroutine eam_splint
1153
1154 subroutine eam_spline(atype, nx, xa, ya, yppa, yp1, ypn, boundary)
1155
1156 ! include 'headers/sizes.h'
1157
1158
1159 ! yp1 and ypn are the first derivatives of y at the two endpoints
1160 ! if boundary is 'L' the lower derivative is used
1161 ! if boundary is 'U' the upper derivative is used
1162 ! if boundary is 'B' then both derivatives are used
1163 ! if boundary is anything else, then both derivatives are assumed to be 0
1164
1165 integer nx, i, k, atype, max_array_size
1166
1167 real( kind = DP ), dimension(:,:) :: xa
1168 real( kind = DP ), dimension(:,:) :: ya
1169 real( kind = DP ), dimension(:,:) :: yppa
1170 real( kind = DP ), allocatable, dimension(:) :: u
1171 real( kind = DP ) :: yp1,ypn,un,qn,sig,p
1172 character boundary
1173
1174 max_array_size = size(xa,1)
1175 allocate(u(max_array_size))
1176
1177
1178 if ((boundary.eq.'l').or.(boundary.eq.'L').or. &
1179 (boundary.eq.'b').or.(boundary.eq.'B')) then
1180 yppa(1, atype) = -0.5E0_DP
1181 u(1) = (3.0E0_DP/(xa(2,atype)-xa(1,atype)))*((ya(2,atype)-&
1182 ya(1,atype))/(xa(2,atype)-xa(1,atype))-yp1)
1183 else
1184 yppa(1,atype) = 0.0E0_DP
1185 u(1) = 0.0E0_DP
1186 endif
1187
1188 do i = 2, nx - 1
1189 sig = (xa(i,atype) - xa(i-1,atype)) / (xa(i+1,atype) - xa(i-1,atype))
1190 p = sig * yppa(i-1,atype) + 2.0E0_DP
1191 yppa(i,atype) = (sig - 1.0E0_DP) / p
1192 u(i) = (6.0E0_DP*((ya(i+1,atype)-ya(i,atype))/(xa(i+1,atype)-xa(i,atype)) - &
1193 (ya(i,atype)-ya(i-1,atype))/(xa(i,atype)-xa(i-1,atype)))/ &
1194 (xa(i+1,atype)-xa(i-1,atype)) - sig * u(i-1))/p
1195 enddo
1196
1197 if ((boundary.eq.'u').or.(boundary.eq.'U').or. &
1198 (boundary.eq.'b').or.(boundary.eq.'B')) then
1199 qn = 0.5E0_DP
1200 un = (3.0E0_DP/(xa(nx,atype)-xa(nx-1,atype)))* &
1201 (ypn-(ya(nx,atype)-ya(nx-1,atype))/(xa(nx,atype)-xa(nx-1,atype)))
1202 else
1203 qn = 0.0E0_DP
1204 un = 0.0E0_DP
1205 endif
1206
1207 yppa(nx,atype)=(un-qn*u(nx-1))/(qn*yppa(nx-1,atype)+1.0E0_DP)
1208
1209 do k = nx-1, 1, -1
1210 yppa(k,atype)=yppa(k,atype)*yppa(k+1,atype)+u(k)
1211 enddo
1212
1213 deallocate(u)
1214 return
1215 end subroutine eam_spline
1216
1217
1218 subroutine calc_eam_rho(r, rho, drho, d2rho, atype)
1219
1220 ! include 'headers/sizes.h'
1221
1222
1223 integer atype, etype, number_r
1224 real( kind = DP ) :: r, rho, drho, d2rho
1225 integer :: i
1226
1227
1228 etype = eam_atype_map(atype)
1229
1230 if (r.lt.eam_rcut(etype)) then
1231 number_r = eam_nr(etype)
1232 call eam_splint(etype, number_r, eam_rvals, eam_rho_r, &
1233 eam_rho_r_pp, r, rho, drho, d2rho)
1234 else
1235 rho = 0.0E0_DP
1236 drho = 0.0E0_DP
1237 d2rho = 0.0E0_DP
1238 endif
1239
1240 return
1241 end subroutine calc_eam_rho
1242
1243 subroutine calc_eam_frho(dens, u, u1, u2, atype)
1244
1245 ! include 'headers/sizes.h'
1246
1247 integer atype, etype, number_rho
1248 real( kind = DP ) :: dens, u, u1, u2
1249 real( kind = DP ) :: rho_vals
1250
1251 etype = eam_atype_map(atype)
1252 number_rho = eam_nrho(etype)
1253 if (dens.lt.eam_rhovals(number_rho, etype)) then
1254 call eam_splint(etype, number_rho, eam_rhovals, eam_f_rho, &
1255 eam_f_rho_pp, dens, u, u1, u2)
1256 else
1257 rho_vals = eam_rhovals(number_rho,etype)
1258 call eam_splint(etype, number_rho, eam_rhovals, eam_f_rho, &
1259 eam_f_rho_pp, rho_vals, u, u1, u2)
1260 endif
1261
1262 return
1263 end subroutine calc_eam_frho
1264
1265 subroutine calc_eam_phi(r, phi, dphi, d2phi, atype)
1266
1267
1268
1269
1270 integer atype, etype, number_r
1271 real( kind = DP ) :: r, phi, dphi, d2phi
1272
1273 etype = eam_atype_map(atype)
1274
1275 if (r.lt.eam_rcut(etype)) then
1276 number_r = eam_nr(etype)
1277 call eam_splint(etype, number_r, eam_rvals, eam_phi_r, &
1278 eam_phi_r_pp, r, phi, dphi, d2phi)
1279 else
1280 phi = 0.0E0_DP
1281 dphi = 0.0E0_DP
1282 d2phi = 0.0E0_DP
1283 endif
1284
1285 return
1286 end subroutine calc_eam_phi
1287
1288
1289 ! Function makes sure atypes are available in potential
1290 function eam_check_type(atype1) result(present)
1291 integer :: atype1
1292 logical :: present
1293 integer :: i, tester
1294
1295 present = .false.
1296
1297 do i = 1, n_eam_atypes
1298 tester = eam_atype(i)
1299 if (atype1 == tester) then
1300 present = .true.
1301 endif
1302 end do
1303
1304 end function eam_check_type
1305
1306
1307
1308 subroutine mass_weight()
1309 integer ia, ja, dim, dim2, idim, idim2
1310 real( kind = DP ) :: mt, m1, m2, wt
1311
1312
1313 do ia = 1, natoms
1314 m1 = mass(ia)
1315 do ja = 1, natoms
1316 m2 = mass(ja)
1317 wt = 1.0E0_DP/dsqrt(m1*m2)
1318 do dim = 1, 3
1319 idim = 3 * (ia-1) + dim
1320 do dim2 = 1, 3
1321 idim2 = 3 * (ja-1) + dim2
1322 d(idim,idim2) = d(idim,idim2)*wt
1323 enddo
1324 enddo
1325 enddo
1326 enddo
1327
1328 end subroutine mass_weight
1329
1330
1331
1332
1333 end module eam_module