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 657 by chuckv, Wed Jul 30 21:17:01 2003 UTC vs.
Revision 843 by mmeineke, Wed Oct 29 20:41:39 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.26 2003-07-30 21:17:01 chuckv Exp $, $Date: 2003-07-30 21:17:01 $, $Name: not supported by cvs2svn $, $Revision: 1.26 $
7 > !! @version $Id: do_Forces.F90,v 1.35 2003-10-29 20:41:39 mmeineke Exp $, $Date: 2003-10-29 20:41:39 $, $Name: not supported by cvs2svn $, $Revision: 1.35 $
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 146 | Line 154 | contains
154  
155  
156      if (FF_uses_EAM) then
157 <       call init_EAM_FF(my_status)
157 >         call init_EAM_FF(my_status)
158         if (my_status /= 0) then
159 +          write(*,*) "init_EAM_FF returned a bad status"
160            thisStat = -1
161            return
162         end if
# Line 206 | Line 215 | contains
215      real( kind = DP ) :: pot_local
216      integer :: nrow
217      integer :: ncol
218 +    integer :: nprocs
219   #endif
220      integer :: nlocal
221      integer :: natoms    
# Line 235 | Line 245 | contains
245      nlocal = getNlocal()
246      natoms = nlocal
247   #endif
248 <    write(*,*) "Inside do_Force Loop"
248 >
249      call check_initialization(localError)
250      if ( localError .ne. 0 ) then
251         call handleError("do_force_loop","Not Initialized")
# Line 247 | Line 257 | contains
257      do_pot = do_pot_c
258      do_stress = do_stress_c
259  
260 +
261      ! Gather all information needed by all force loops:
262      
263   #ifdef IS_MPI    
# Line 263 | Line 274 | contains
274      endif
275      
276   #endif
277 <    
277 >
278 > !! Begin force loop timing:
279 > #ifdef PROFILE
280 >    call cpu_time(forceTimeInitial)
281 >    nloops = nloops + 1
282 > #endif
283 >  
284      if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then
285         !! See if we need to update neighbor lists
286         call checkNeighborList(nlocal, q, listSkin, update_nlist)  
# Line 271 | Line 288 | contains
288         !! do_prepair_loop_if_needed
289         !! if_mpi_scatter_stuff_from_prepair
290         !! if_mpi_gather_stuff_from_prepair_to_main_loop
291 <
291 >    
292   !--------------------PREFORCE LOOP----------->>>>>>>>>>>>>>>>>>>>>>>>>>>
293   #ifdef IS_MPI
294      
# Line 347 | Line 364 | contains
364         neighborListSize = size(list)
365    
366         nlist = 0
367 <      
367 >
368         do i = 1, natoms-1
369            point(i) = nlist + 1
370            
# Line 359 | Line 376 | contains
376            
377  
378               if (rijsq < rlistsq) then
379 <                
379 >
380 >          
381                  nlist = nlist + 1
382                
383                  if (nlist > neighborListSize) then
# Line 384 | Line 402 | contains
402         point(natoms) = nlist + 1
403        
404      else !! (update)
405 <      
405 >  
406         ! use the list to find the neighbors
407         do i = 1, natoms-1
408            JBEG = POINT(i)
# Line 405 | Line 423 | contains
423      endif    
424   #endif
425      !! Do rest of preforce calculations
426 <   call do_preforce(nlocal,pot)
426 >    !! do necessary preforce calculations  
427 >    call do_preforce(nlocal,pot)
428 >   ! we have already updated the neighbor list set it to false...
429 >   update_nlist = .false.
430      else
431         !! See if we need to update neighbor lists for non pre-pair
432         call checkNeighborList(nlocal, q, listSkin, update_nlist)  
# Line 558 | Line 579 | contains
579   #endif
580      
581      ! phew, done with main loop.
582 <    
582 >
583 > !! Do timing
584 > #ifdef PROFILE
585 >    call cpu_time(forceTimeFinal)
586 >    forceTime = forceTime + forceTimeFinal - forceTimeInitial
587 > #endif
588 >
589 >
590   #ifdef IS_MPI
591      !!distribute forces
592    
# Line 672 | Line 700 | contains
700      endif
701  
702   #endif
703 <    
703 >
704 > #ifdef PROFILE
705 >    if (do_pot) then
706 >
707 > #ifdef IS_MPI
708 >
709 >      
710 >       call printCommTime()
711 >
712 >       call mpi_allreduce(forceTime,globalForceTime,1,MPI_DOUBLE_PRECISION, &
713 >            mpi_sum,mpi_comm_world,mpi_err)
714 >
715 >       call mpi_allreduce(forceTime,maxForceTime,1,MPI_DOUBLE_PRECISION, &
716 >            MPI_MAX,mpi_comm_world,mpi_err)
717 >      
718 >       call mpi_comm_size( MPI_COMM_WORLD, nprocs,mpi_err)
719 >      
720 >       if (getMyNode() == 0) then
721 >          write(*,*) "Total processor time spent in force calculations is: ", globalForceTime
722 >          write(*,*) "Total Time spent in force loop per processor is: ", globalforceTime/nprocs
723 >          write(*,*) "Maximum force time on any processor is: ", maxForceTime
724 >       end if
725 > #else
726 >       write(*,*) "Time spent in force loop is: ", forceTime
727 > #endif
728 >
729 >    
730 >    endif
731 >
732 > #endif
733 >
734    end subroutine do_force_loop
735  
736    subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot)
# Line 725 | Line 783 | contains
783         call getElementProperty(atypes, me_j, "is_DP", is_DP_j)
784        
785         if ( is_DP_i .and. is_DP_j ) then
728          
786            call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, &
787                 do_pot, do_stress)
788            if (FF_uses_RF .and. SimUsesRF()) then
# Line 750 | Line 807 | contains
807  
808      if (FF_uses_GB .and. SimUsesGB()) then
809  
810 +
811         call getElementProperty(atypes, me_i, "is_GB", is_GB_i)
812         call getElementProperty(atypes, me_j, "is_GB", is_GB_j)
813        
# Line 795 | Line 853 | contains
853    
854     r = sqrt(rijsq)
855    
856 +
857   #ifdef IS_MPI
858     if (tagRow(i) .eq. tagColumn(j)) then
859        write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
# Line 809 | Line 868 | contains
868     me_j = atid(j)
869    
870   #endif
871 <  
871 >    
872     if (FF_uses_EAM .and. SimUsesEAM()) then
873        call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
874        call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
# Line 817 | Line 876 | contains
876        if ( is_EAM_i .and. is_EAM_j ) &
877             call calc_EAM_prepair_rho(i, j, d, r, rijsq )
878     endif
820  end subroutine do_prepair
879  
880 + end subroutine do_prepair
881  
882  
883  
884 +
885    subroutine do_preforce(nlocal,pot)
886      integer :: nlocal
887      real( kind = dp ) :: pot
888  
889 <   if (FF_uses_EAM .and. SimUsesEAM()) then
890 <      call calc_EAM_preforce_Frho(nlocal,pot)
891 <   endif
889 >    if (FF_uses_EAM .and. SimUsesEAM()) then
890 >       call calc_EAM_preforce_Frho(nlocal,pot)
891 >    endif
892  
893  
894    end subroutine do_preforce
# Line 939 | Line 999 | contains
999  
1000   #endif
1001  
1002 +
1003 +    if (FF_uses_EAM .and. SimUsesEAM()) then
1004 +       call clean_EAM()
1005 +    endif
1006 +
1007 +
1008 +
1009 +
1010 +
1011      rf = 0.0_dp
1012      tau_Temp = 0.0_dp
1013      virial_Temp = 0.0_dp
# Line 1051 | Line 1120 | end module do_Forces
1120      doesit = FF_uses_RF
1121    end function FF_RequiresPostpairCalc
1122    
1123 + !! This cleans componets of force arrays belonging only to fortran
1124 +
1125   end module do_Forces

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines