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 619 by mmeineke, Tue Jul 15 22:22:41 2003 UTC vs.
Revision 626 by mmeineke, Wed Jul 16 21:30:56 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.20 2003-07-15 22:22:41 mmeineke Exp $, $Date: 2003-07-15 22:22:41 $, $Name: not supported by cvs2svn $, $Revision: 1.20 $
7 > !! @version $Id: do_Forces.F90,v 1.21 2003-07-16 21:30:55 mmeineke Exp $, $Date: 2003-07-16 21:30:55 $, $Name: not supported by cvs2svn $, $Revision: 1.21 $
8  
9   module do_Forces
10    use force_globals
# Line 17 | Line 17 | module do_Forces
17    use dipole_dipole
18    use reaction_field
19    use gb_pair
20 +  use vector_class
21   #ifdef IS_MPI
22    use mpiSimulation
23   #endif
# Line 27 | Line 28 | module do_Forces
28   #define __FORTRAN90
29   #include "fForceField.h"
30  
31 <  logical, save :: do_forces_initialized = .false.
31 >  logical, save :: do_forces_initialized = .false., haveRlist = .false.
32 >  logical, save :: havePolicies = .false.
33    logical, save :: FF_uses_LJ
34    logical, save :: FF_uses_sticky
35    logical, save :: FF_uses_dipoles
# Line 35 | Line 37 | module do_Forces
37    logical, save :: FF_uses_GB
38    logical, save :: FF_uses_EAM
39  
40 +  real(kind=dp), save :: rlist, rlistsq
41 +
42    public :: init_FF
43    public :: do_force_loop
44 +  public :: setRlistDF
45  
46   contains
47  
48 +  subroutine setRlistDF( this_rlist )
49 +    
50 +    real(kind=dp) :: this_rlist
51 +
52 +    rlist = this_rlist
53 +    rlistsq = rlist * rlist
54 +    
55 +    haveRlist = .true.
56 +    if( havePolicies ) do_forces_initialized = .true.
57 +
58 +  end subroutine setRlistDF    
59 +
60    subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat)
61  
62      integer, intent(in) :: LJMIXPOLICY
# Line 87 | Line 104 | contains
104      !! check to make sure the FF_uses_RF setting makes sense
105      
106      if (FF_uses_dipoles) then
90       rrf = getRrf()
91       rt = getRt()      
92       call initialize_dipole(rrf, rt)
107         if (FF_uses_RF) then
108            dielect = getDielect()
109 <          call initialize_rf(rrf, rt, dielect)
109 >          call initialize_rf(dielect)
110         endif
111      else
112         if (FF_uses_RF) then          
# Line 100 | Line 114 | contains
114            thisStat = -1
115            return
116         endif
117 <    endif
117 >    endif
118  
119      if (FF_uses_LJ) then
120        
107       call getRcut(rcut)
108
121         select case (LJMIXPOLICY)
122         case (LB_MIXING_RULE)
123 <          call init_lj_FF(LB_MIXING_RULE, rcut, my_status)            
123 >          call init_lj_FF(LB_MIXING_RULE, my_status)            
124         case (EXPLICIT_MIXING_RULE)
125 <          call init_lj_FF(EXPLICIT_MIXING_RULE, rcut, my_status)
125 >          call init_lj_FF(EXPLICIT_MIXING_RULE, my_status)
126         case default
127            write(default_error,*) 'unknown LJ Mixing Policy!'
128            thisStat = -1
# Line 150 | Line 162 | contains
162         endif
163      endif
164  
165 <    do_forces_initialized = .true.    
166 <
165 >    havePolicies = .true.
166 >    if( haveRlist ) do_forces_initialized = .true.
167 >    
168    end subroutine init_FF
169    
170  
# Line 185 | Line 198 | contains
198      logical :: update_nlist  
199      integer :: i, j, jbeg, jend, jnab
200      integer :: nlist
201 <    real( kind = DP ) ::  rijsq, rlistsq, rcutsq, rlist, rcut, rrf, rt, dielect
201 >    real( kind = DP ) ::  rijsq
202      real(kind=dp),dimension(3) :: d
203      real(kind=dp) :: rfpot, mu_i, virial
204      integer :: me_i
# Line 193 | Line 206 | contains
206      integer :: neighborListSize
207      integer :: listerror, error
208      integer :: localError
209 +
210 +    real(kind=dp) :: listSkin = 1.0
211      
212  
213      !! initialize local variables  
# Line 207 | Line 222 | contains
222      natoms = nlocal
223   #endif
224    
210    call getRcut(rcut,rc2=rcutsq)
211    call getRlist(rlist,rlistsq)
212    rt = getRt()
213    rrf = getRrf()
214    dielect = getDielect()
215    
216    if( FF_uses_LJ) then      
217       call lj_new_rcut( rcut, localError )
218       if ( localError .ne. 0 ) then
219          error = -1
220          return
221       end if
222    end if
223    
224    
225    if( FF_uses_dipoles ) then
226      
227       if( rcut .lt. rrf ) then
228          rcut = rrf
229          rlist = rcut + 1.0_dp
230          rcutsq = rcut * rcut
231          rlistsq = rlist * rlist
232       end if
233      
234       call initialize_dipole( rrf, rt )
235    end if
236    
237    if( FF_uses_RF )call initialize_rf( rrf, rt, dielect )
238    
239    
225      call check_initialization(localError)
226      if ( localError .ne. 0 ) then
227         error = -1
# Line 266 | Line 251 | contains
251      
252      if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then
253         !! See if we need to update neighbor lists
254 <       call checkNeighborList(nlocal, q, rcut, rlist, update_nlist)  
254 >       call checkNeighborList(nlocal, q, listSkin, update_nlist)  
255         !! if_mpi_gather_stuff_for_prepair
256         !! do_prepair_loop_if_needed
257         !! if_mpi_scatter_stuff_from_prepair
258         !! if_mpi_gather_stuff_from_prepair_to_main_loop
259      else
260         !! See if we need to update neighbor lists
261 <       call checkNeighborList(nlocal, q, rcut, rlist, update_nlist)  
261 >       call checkNeighborList(nlocal, q, listSkin, update_nlist)  
262      endif
263      
264   #ifdef IS_MPI
# Line 296 | Line 281 | contains
281              
282               call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq)
283              
284 <             if (rijsq <  rlistsq) then            
284 >             if (rijsq < rlistsq) then            
285                  
286                  nlist = nlist + 1
287                  
# Line 312 | Line 297 | contains
297                  
298                  list(nlist) = j
299                                  
300 <                if (rijsq <  rcutsq) then
301 <                   call do_pair(i, j, rijsq, d, do_pot, do_stress, &
302 <                        u_l, A, f, t, pot_local)
318 <                endif
300 >                call do_pair(i, j, rijsq, d, do_pot, do_stress, &
301 >                     u_l, A, f, t, pot_local)
302 >                
303               endif
304            enddo inner
305         enddo
# Line 365 | Line 349 | contains
349               call get_interatomic_vector(q(:,i), q(:,j), d, rijsq)
350            
351  
352 <             if (rijsq <  rlistsq) then
352 >             if (rijsq < rlistsq) then
353                  
354                  nlist = nlist + 1
355                
# Line 381 | Line 365 | contains
365                  
366                  list(nlist) = j
367                  
368 <                if (rijsq <  rcutsq) then
385 <                   call do_pair(i, j, rijsq, d, do_pot, do_stress, &
368 >                call do_pair(i, j, rijsq, d, do_pot, do_stress, &
369                          u_l, A, f, t, pot)
370 <                endif
370 >                
371               endif
372            enddo inner
373         enddo

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines