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 845 by gezelter, Thu Oct 30 18:59:20 2003 UTC vs.
Revision 883 by chuckv, Thu Dec 18 20:46:45 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.36 2003-10-30 18:59:20 gezelter Exp $, $Date: 2003-10-30 18:59:20 $, $Name: not supported by cvs2svn $, $Revision: 1.36 $
7 > !! @version $Id: do_Forces.F90,v 1.38 2003-12-18 20:46:45 chuckv Exp $, $Date: 2003-12-18 20:46:45 $, $Name: not supported by cvs2svn $, $Revision: 1.38 $
8  
9   module do_Forces
10    use force_globals
# Line 46 | Line 46 | module do_Forces
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
49 >  public :: getforcetime
50 >  real, save :: forceTime = 0
51 >  real :: forceTimeInitial, forceTimeFinal
52   #endif
53  
54   contains
# Line 444 | Line 442 | contains
442   #ifdef IS_MPI
443      
444      if (update_nlist) then
447      
445         !! save current configuration, construct neighbor list,
446         !! and calculate forces
447         call saveNeighborList(nlocal, q)
# Line 510 | Line 507 | contains
507   #else
508      
509      if (update_nlist) then
510 <      
510 >
511         ! save current configuration, contruct neighbor list,
512         ! and calculate forces
513         call saveNeighborList(natoms, q)
# Line 700 | Line 697 | contains
697  
698   #endif
699  
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
700      
701      endif
702  
# Line 1119 | Line 1091 | contains
1091      doesit = FF_uses_RF
1092    end function FF_RequiresPostpairCalc
1093    
1094 + #ifdef PROFILE
1095 +  function getforcetime() return(totalforcetime)
1096 +    real(kind=dp) :: totalforcetime
1097 +    totalforcetime = forcetime
1098 +  end function getforcetime
1099 + #endif
1100 +
1101   !! This cleans componets of force arrays belonging only to fortran
1102  
1103   end module do_Forces

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines