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 673 by chuckv, Fri Aug 8 21:22:37 2003 UTC vs.
Revision 694 by chuckv, Wed Aug 13 21:20:20 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.28 2003-08-08 21:22:37 chuckv Exp $, $Date: 2003-08-08 21:22:37 $, $Name: not supported by cvs2svn $, $Revision: 1.28 $
7 > !! @version $Id: do_Forces.F90,v 1.29 2003-08-13 21:20:20 chuckv Exp $, $Date: 2003-08-13 21:20:20 $, $Name: not supported by cvs2svn $, $Revision: 1.29 $
8  
9   module do_Forces
10    use force_globals
# Line 45 | Line 45 | contains
45    public :: do_force_loop
46    public :: setRlistDF
47  
48 + #ifdef PROFILE
49 +  real(kind = dp) :: forceTime
50 +  real(kind = dp) :: forceTimeInitial, forceTimeFinal
51 +  real(kind = dp) :: globalForceTime
52 +  real(kind = dp) :: maxForceTime
53 +  integer, save :: nloops = 0
54 + #endif
55 +
56   contains
57  
58    subroutine setRlistDF( this_rlist )
# Line 206 | Line 214 | contains
214      real( kind = DP ) :: pot_local
215      integer :: nrow
216      integer :: ncol
217 +    integer :: nprocs
218   #endif
219      integer :: nlocal
220      integer :: natoms    
# Line 263 | Line 272 | contains
272         call gather(A,A_Col,plan_col_rotation)
273      endif
274      
275 + #endif
276 +
277 + !! Begin force loop timing:
278 + #ifdef PROFILE
279 +    call cpu_time(forceTimeInitial)
280 +    nloops = nloops + 1
281   #endif
282    
283      if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then
# Line 563 | Line 578 | contains
578   #endif
579      
580      ! phew, done with main loop.
581 <    
581 >
582 > !! Do timing
583 > #ifdef PROFILE
584 >    call cpu_time(forceTimeFinal)
585 >    forceTime = forceTime + forceTimeFinal - forceTimeInitial
586 > #endif
587 >
588 >
589   #ifdef IS_MPI
590      !!distribute forces
591    
# Line 677 | Line 699 | contains
699      endif
700  
701   #endif
702 <    
702 >
703 > #ifdef PROFILE
704 >    if (do_pot) then
705 >
706 > #ifdef IS_MPI
707 >
708 >      
709 >       call printCommTime()
710 >
711 >       call mpi_allreduce(forceTime,globalForceTime,1,MPI_DOUBLE_PRECISION, &
712 >            mpi_sum,mpi_comm_world,mpi_err)
713 >
714 >       call mpi_allreduce(forceTime,maxForceTime,1,MPI_DOUBLE_PRECISION, &
715 >            MPI_MAX,mpi_comm_world,mpi_err)
716 >      
717 >       call mpi_comm_size( MPI_COMM_WORLD, nprocs,mpi_err)
718 >      
719 >       if (getMyNode() == 0) then
720 >          write(*,*) "Total processor time spent in force calculations is: ", globalForceTime
721 >          write(*,*) "Total Time spent in force loop per processor is: ", globalforceTime/nprocs
722 >          write(*,*) "Maximum force time on any processor is: ", maxForceTime
723 >       end if
724 > #else
725 >       write(*,*) "Time spent in force loop is: ", forceTime
726 > #endif
727 >
728 >    
729 >    endif
730 >
731 > #endif
732 >
733 >
734 >
735    end subroutine do_force_loop
736  
737    subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines