ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_eam.F90
(Generate patch)

Comparing trunk/OOPSE/libmdtools/calc_eam.F90 (file contents):
Revision 801 by chuckv, Sat Oct 4 18:46:12 2003 UTC vs.
Revision 882 by chuckv, Wed Dec 17 20:13:33 2003 UTC

# Line 66 | Line 66 | module eam
66    real( kind = dp),save, dimension(:), allocatable :: frho_col
67    real( kind = dp),save, dimension(:), allocatable :: rho_row
68    real( kind = dp),save, dimension(:), allocatable :: rho_col
69 +  real( kind = dp),save, dimension(:), allocatable :: rho_tmp
70    real( kind = dp),save, dimension(:), allocatable :: d2frhodrhodrho_col
71    real( kind = dp),save, dimension(:), allocatable :: d2frhodrhodrho_row
72   #endif
# Line 222 | Line 223 | contains
223                 0.0E0_DP, 0.0E0_DP, 'N')
224         enddo
225  
225
226   !       current_rcut_max = EAMList%EAMParams(1)%eam_rcut
227         !! find the smallest rcut for any eam atype
228   !       do i = 2, EAMList%currentAddition
# Line 294 | Line 294 | contains
294      
295   #ifdef IS_MPI
296  
297 +    if (allocated(rho_tmp)) deallocate(rho_tmp)
298 +    allocate(rho_tmp(nlocal),stat=alloc_stat)
299 +    if (alloc_stat /= 0) then
300 +       status = -1
301 +       return
302 +    end if
303 +
304 +
305      if (allocated(frho_row)) deallocate(frho_row)
306      allocate(frho_row(nrow),stat=alloc_stat)
307      if (alloc_stat /= 0) then
# Line 377 | Line 385 | contains
385      frho_col = 0.0_dp
386      rho_row  = 0.0_dp
387      rho_col  = 0.0_dp
388 +    rho_tmp  = 0.0_dp
389      dfrhodrho_row = 0.0_dp
390      dfrhodrho_col = 0.0_dp
391   #endif
# Line 525 | Line 534 | contains
534                 EAMList%EAMParams(myid_atom2)%eam_rho_r_pp, &
535                 r, rho_j_at_i,drho,d2rho)
536  
537 <
537 >    
538        
539        
540   #ifdef  IS_MPI
# Line 536 | Line 545 | contains
545         endif
546  
547  
548 +
549 +
550 +
551 +
552    end subroutine calc_eam_prepair_rho
553  
554  
# Line 561 | Line 574 | contains
574        write(errMsg,*) " Error scattering rho_row into rho"
575        call handleError(RoutineName,errMesg)
576     endif      
577 <    call scatter(rho_col,rho,plan_col,eam_err)
577 >    call scatter(rho_col,rho_tmp,plan_col,eam_err)
578      if (eam_err /= 0 ) then
579        write(errMsg,*) " Error scattering rho_col into rho"
580        call handleError(RoutineName,errMesg)
581     endif
582 +
583 +      rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
584   #endif
585  
586  
587 +
588   !! Calculate F(rho) and derivative
589      do atom = 1, nlocal
590         me = atid(atom)

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines