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 669 by chuckv, Thu Aug 7 00:47:33 2003 UTC vs.
Revision 673 by chuckv, Fri Aug 8 21:22:37 2003 UTC

# Line 51 | Line 51 | module eam
51  
52  
53    !! Arrays for derivatives used in force calculation
54 <  real( kind = dp),save, dimension(:), allocatable :: frho
55 <  real( kind = dp),save, dimension(:), allocatable :: rho
54 >  real( kind = dp), dimension(:), allocatable :: frho
55 >  real( kind = dp), dimension(:), allocatable :: rho
56  
57 <  real( kind = dp),save, dimension(:), allocatable :: dfrhodrho
58 <  real( kind = dp),save, dimension(:), allocatable :: d2frhodrhodrho
57 >  real( kind = dp), dimension(:), allocatable :: dfrhodrho
58 >  real( kind = dp), dimension(:), allocatable :: d2frhodrhodrho
59  
60  
61   !! Arrays for MPI storage
# Line 89 | Line 89 | module eam
89    public :: newEAMtype
90    public :: calc_eam_prepair_rho
91    public :: calc_eam_preforce_Frho
92 <  
92 >  public :: clean_EAM
93  
94   contains
95  
# Line 233 | Line 233 | contains
233   !       do i = 1, EAMList%currentAddition
234   !          EAMList%EAMParam(i)s%eam_atype_map(eam_atype(i)) = i
235   !       end do
236
237
236         !! Allocate arrays for force calculation
237            call allocateEAM(alloc_stat)
238            if (alloc_stat /= 0 ) then
# Line 365 | Line 363 | contains
363  
364  
365    subroutine clean_EAM()
366 +  
367     ! clean non-IS_MPI first
368      frho = 0.0_dp
369      rho  = 0.0_dp
# Line 484 | Line 483 | contains
483      integer :: myid_atom2
484  
485   ! check to see if we need to be cleaned at the start of a force loop
487    
488    if (cleanme) then
489       call clean_EAM
490       cleanme = .false.
491    end if
486        
487  
494    
488  
489 +
490   #ifdef IS_MPI
491      myid_atom1 = atid_Row(atom1)
492      myid_atom2 = atid_Col(atom2)
# Line 502 | Line 496 | contains
496   #endif
497  
498      if (r.lt.EAMList%EAMParams(myid_atom1)%eam_rcut) then
505
499  
500  
501 +
502         call eam_splint(EAMList%EAMParams(myid_atom1)%eam_nr, &
503              EAMList%EAMParams(myid_atom1)%eam_rvals, &
504              EAMList%EAMParams(myid_atom1)%eam_rho_r, &
# Line 518 | Line 512 | contains
512   #else
513         rho(atom2) = rho(atom2) + rho_i_at_j
514   #endif
515 + !       write(*,*) atom1,atom2,r,rho_i_at_j
516         endif
517  
518         if (r.lt.EAMList%EAMParams(myid_atom2)%eam_rcut) then
# Line 537 | Line 532 | contains
532   #endif
533         endif
534  
540    
535  
536    end subroutine calc_eam_prepair_rho
537  
# Line 554 | Line 548 | contains
548      integer :: atype1
549      integer :: me
550      integer :: n_rho_points
551 <    ! reset clean forces to be true at top of calc rho.
558 <    cleanme = .true.
551 >
552    
553 < !! Scatter the electron density from  pre-pair calculation back to local atoms
553 >    cleanme = .true.
554 >    !! Scatter the electron density from  pre-pair calculation back to local atoms
555   #ifdef IS_MPI
556      call scatter(rho_row,rho,plan_row,eam_err)
557      if (eam_err /= 0 ) then
# Line 595 | Line 589 | contains
589         end if
590  
591  
592 <       frho(i) = u
593 <       dfrhodrho(i) = u1
594 <       d2frhodrhodrho(i) = u2
592 >       frho(atom) = u
593 >       dfrhodrho(atom) = u1
594 >       d2frhodrhodrho(atom) = u2
595         pot = pot + u
596 +
597      enddo
598  
599 <  
599 >  
600  
601   #ifdef IS_MPI
602      !! communicate f(rho) and derivatives back into row and column arrays
# Line 632 | Line 627 | contains
627      endif
628   #endif
629  
630 +  
631    end subroutine calc_eam_preforce_Frho
632  
633  
# Line 663 | Line 659 | contains
659      integer :: id1,id2
660      integer  :: mytype_atom1
661      integer  :: mytype_atom2
666
662  
663   !Local Variables
664      
665 <
665 >   ! write(*,*) "Frho: ", Frho(atom1)
666  
667      phab = 0.0E0_DP
668      dvpdr = 0.0E0_DP

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines