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

Comparing trunk/OOPSE/libmdtools/do_Forces.F90 (file contents):
Revision 441 by chuckv, Tue Apr 1 16:50:14 2003 UTC vs.
Revision 480 by chuckv, Tue Apr 8 17:16:22 2003 UTC

# Line 4 | Line 4
4  
5   !! @author Charles F. Vardeman II
6   !! @author Matthew Meineke
7 < !! @version $Id: do_Forces.F90,v 1.7 2003-04-01 16:50:14 chuckv Exp $, $Date: 2003-04-01 16:50:14 $, $Name: not supported by cvs2svn $, $Revision: 1.7 $
7 > !! @version $Id: do_Forces.F90,v 1.12 2003-04-08 17:16:22 chuckv Exp $, $Date: 2003-04-08 17:16:22 $, $Name: not supported by cvs2svn $, $Revision: 1.12 $
8  
9   module do_Forces
10    use force_globals
# Line 140 | Line 140 | contains
140  
141      if (FF_uses_GB .and. FF_uses_LJ) then
142      endif
143 <
143 >    if (.not. do_forces_initialized) then
144 >       !! Create neighbor lists
145 >       call expandNeighborList(getNlocal(), my_status)
146 >       if (my_Status /= 0) then
147 >          write(default_error,*) "SimSetup: ExpandNeighborList returned error."
148 >          thisStat = -1
149 >          return
150 >       endif
151 >    endif
152  
153      do_forces_initialized = .true.    
154  
# Line 210 | Line 218 | contains
218  
219      do_pot = do_pot_c
220      do_stress = do_stress_c
221 +    
222  
223      ! Gather all information needed by all force loops:
224      
# Line 246 | Line 255 | contains
255        
256         !! save current configuration, construct neighbor list,
257         !! and calculate forces
258 <       call saveNeighborList(q)
258 >       call saveNeighborList(nlocal, q)
259        
260         neighborListSize = size(list)
261         nlist = 0      
# Line 313 | Line 322 | contains
322        
323         ! save current configuration, contruct neighbor list,
324         ! and calculate forces
325 <       call saveNeighborList(q)
325 >       call saveNeighborList(natoms, q)
326        
327         neighborListSize = size(list)
328    
# Line 479 | Line 488 | contains
488      endif
489  
490      if (do_stress) then
491 <       call mpi_allreduce(tau, tau_Temp,9,mpi_double_precision,mpi_sum, &
491 >       call mpi_allreduce(tau_Temp, tau,9,mpi_double_precision,mpi_sum, &
492              mpi_comm_world,mpi_err)
493 <       call mpi_allreduce(virial, virial_Temp,1,mpi_double_precision,mpi_sum, &
493 >       call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
494              mpi_comm_world,mpi_err)
495      endif
496  
# Line 499 | Line 508 | contains
508    subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
509  
510      real( kind = dp ) :: pot
511 <    real( kind = dp ), dimension(:,:) :: u_l
512 <    real (kind=dp), dimension(:,:) :: A
513 <    real (kind=dp), dimension(:,:) :: f
514 <    real (kind=dp), dimension(:,:) :: t
511 >    real( kind = dp ), dimension(3,getNlocal()) :: u_l
512 >    real (kind=dp), dimension(9,getNlocal()) :: A
513 >    real (kind=dp), dimension(3,getNlocal()) :: f
514 >    real (kind=dp), dimension(3,getNlocal()) :: t
515  
516      logical, intent(inout) :: do_pot, do_stress
517      integer, intent(in) :: i, j
# Line 655 | Line 664 | contains
664      rf = 0.0_dp
665      tau_Temp = 0.0_dp
666      virial_Temp = 0.0_dp
658    
667    end subroutine zero_work_arrays
668    
669    function skipThisPair(atom1, atom2) result(skip_it)

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines