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 648 by chuckv, Wed Jul 23 22:13:59 2003 UTC vs.
Revision 673 by chuckv, Fri Aug 8 21:22:37 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.24 2003-07-23 22:13:59 chuckv Exp $, $Date: 2003-07-23 22:13:59 $, $Name: not supported by cvs2svn $, $Revision: 1.24 $
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 $
8  
9   module do_Forces
10    use force_globals
# Line 18 | Line 18 | module do_Forces
18    use reaction_field
19    use gb_pair
20    use vector_class
21 +  use eam
22 +  use status
23   #ifdef IS_MPI
24    use mpiSimulation
25   #endif
# Line 141 | Line 143 | contains
143            return
144         end if
145      endif
146 +
147 +
148 +    if (FF_uses_EAM) then
149 +       call init_EAM_FF(my_status)
150 +       if (my_status /= 0) then
151 +          thisStat = -1
152 +          return
153 +       end if
154 +    endif
155 +
156 +
157      
158      if (FF_uses_GB) then
159         call check_gb_pair_FF(my_status)
# Line 161 | Line 174 | contains
174            return
175         endif
176      endif
177 +    
178  
179      havePolicies = .true.
180      if( haveRlist ) do_forces_initialized = .true.
181 <    
181 >
182    end subroutine init_FF
183    
184  
# Line 221 | Line 235 | contains
235      nlocal = getNlocal()
236      natoms = nlocal
237   #endif
238 <  
238 >
239      call check_initialization(localError)
240      if ( localError .ne. 0 ) then
241 +       call handleError("do_force_loop","Not Initialized")
242         error = -1
243         return
244      end if
# Line 232 | Line 247 | contains
247      do_pot = do_pot_c
248      do_stress = do_stress_c
249  
250 +
251      ! Gather all information needed by all force loops:
252      
253   #ifdef IS_MPI    
# Line 248 | Line 264 | contains
264      endif
265      
266   #endif
267 <    
267 >  
268      if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then
269         !! See if we need to update neighbor lists
270         call checkNeighborList(nlocal, q, listSkin, update_nlist)  
# Line 256 | Line 272 | contains
272         !! do_prepair_loop_if_needed
273         !! if_mpi_scatter_stuff_from_prepair
274         !! if_mpi_gather_stuff_from_prepair_to_main_loop
275 <
275 >    
276   !--------------------PREFORCE LOOP----------->>>>>>>>>>>>>>>>>>>>>>>>>>>
277   #ifdef IS_MPI
278      
# Line 332 | Line 348 | contains
348         neighborListSize = size(list)
349    
350         nlist = 0
351 <      
351 >
352         do i = 1, natoms-1
353            point(i) = nlist + 1
354            
# Line 344 | Line 360 | contains
360            
361  
362               if (rijsq < rlistsq) then
363 <                
363 >
364 >          
365                  nlist = nlist + 1
366                
367                  if (nlist > neighborListSize) then
# Line 369 | Line 386 | contains
386         point(natoms) = nlist + 1
387        
388      else !! (update)
389 <      
389 >  
390         ! use the list to find the neighbors
391         do i = 1, natoms-1
392            JBEG = POINT(i)
# Line 390 | Line 407 | contains
407      endif    
408   #endif
409      !! Do rest of preforce calculations
410 <   call do_preforce(nlocal,pot)
410 >    !! do necessary preforce calculations  
411 >    call do_preforce(nlocal,pot)
412 >   ! we have already updated the neighbor list set it to false...
413 >   update_nlist = .false.
414      else
415         !! See if we need to update neighbor lists for non pre-pair
416         call checkNeighborList(nlocal, q, listSkin, update_nlist)  
# Line 780 | Line 800 | contains
800    
801     r = sqrt(rijsq)
802    
803 +
804   #ifdef IS_MPI
805     if (tagRow(i) .eq. tagColumn(j)) then
806        write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
# Line 794 | Line 815 | contains
815     me_j = atid(j)
816    
817   #endif
818 <  
818 >    
819     if (FF_uses_EAM .and. SimUsesEAM()) then
820        call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i)
821        call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j)
# Line 802 | Line 823 | contains
823        if ( is_EAM_i .and. is_EAM_j ) &
824             call calc_EAM_prepair_rho(i, j, d, r, rijsq )
825     endif
805  end subroutine do_prepair
826  
827 + end subroutine do_prepair
828  
829  
830  
831 +
832    subroutine do_preforce(nlocal,pot)
833      integer :: nlocal
834      real( kind = dp ) :: pot
835  
836 <   if (FF_uses_EAM .and. SimUsesEAM()) then
837 <      call calc_EAM_preforce_Frho(nlocal,pot)
838 <   endif
836 >    if (FF_uses_EAM .and. SimUsesEAM()) then
837 >       call calc_EAM_preforce_Frho(nlocal,pot)
838 >    endif
839  
840  
841    end subroutine do_preforce
# Line 876 | Line 898 | contains
898      error = 0
899      ! Make sure we are properly initialized.
900      if (.not. do_forces_initialized) then
901 +       write(*,*) "Forces not initialized"
902         error = -1
903         return
904      endif
# Line 923 | Line 946 | contains
946  
947   #endif
948  
949 +
950 +    if (FF_uses_EAM .and. SimUsesEAM()) then
951 +       call clean_EAM()
952 +    endif
953 +
954 +
955 +
956 +
957 +
958      rf = 0.0_dp
959      tau_Temp = 0.0_dp
960      virial_Temp = 0.0_dp
# Line 1035 | Line 1067 | end module do_Forces
1067      doesit = FF_uses_RF
1068    end function FF_RequiresPostpairCalc
1069    
1070 + !! This cleans componets of force arrays belonging only to fortran
1071 +
1072   end module do_Forces

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines