ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_LJ_FF.F90
(Generate patch)

Comparing trunk/OOPSE/libmdtools/calc_LJ_FF.F90 (file contents):
Revision 483 by gezelter, Wed Apr 9 04:06:43 2003 UTC vs.
Revision 1150 by gezelter, Fri May 7 21:35:05 2004 UTC

# Line 2 | Line 2
2   !! Corresponds to the force field defined in lj_FF.cpp
3   !! @author Charles F. Vardeman II
4   !! @author Matthew Meineke
5 < !! @version $Id: calc_LJ_FF.F90,v 1.4 2003-04-09 04:06:43 gezelter Exp $, $Date: 2003-04-09 04:06:43 $, $Name: not supported by cvs2svn $, $Revision: 1.4 $
5 > !! @version $Id: calc_LJ_FF.F90,v 1.18 2004-05-07 21:35:04 gezelter Exp $, $Date: 2004-05-07 21:35:04 $, $Name: not supported by cvs2svn $, $Revision: 1.18 $
6  
7   module lj
8    use definitions
9    use atype_module
10 +  use switcheroo
11    use vector_class
12    use simulation
13   #ifdef IS_MPI
# Line 21 | Line 22 | module lj
22   #include "fForceField.h"
23  
24    integer, save :: LJ_Mixing_Policy
25 <  integer, save :: LJ_rcut
25 >  real(kind=DP), save :: LJ_rcut
26 >  logical, save :: havePolicy = .false., haveCut = .false.
27 >
28    
29    !! Logical has lj force field module been initialized?
30    
# Line 29 | Line 32 | module lj
32    
33    !! Public methods and data
34    public :: init_LJ_FF
35 <  public :: LJ_new_rcut
35 >  public :: setCutoffLJ
36    public :: do_lj_pair
37    
38    type :: lj_mixed_params
# Line 45 | Line 48 | module lj
48  
49       real ( kind = dp )  :: delta  = 0.0_dp
50  
51 +
52    end type lj_mixed_params
53    
54    type (lj_mixed_params), dimension(:,:), pointer :: ljMixed
55    
56   contains
57  
58 <  subroutine init_LJ_FF(mix_Policy, rcut, status)
58 >  subroutine init_LJ_FF(mix_Policy, status)
59      integer, intent(in) :: mix_Policy
60      integer, intent(out) :: status
61      integer :: myStatus
58    real(kind=dp) :: rcut
62      
63      if (mix_Policy == LB_MIXING_RULE) then
64         LJ_Mixing_Policy = LB_MIXING_RULE
# Line 68 | Line 71 | contains
71            return
72         endif
73      endif
71    
72    LJ_rcut = rcut
74  
75 <    status = 0
76 <    call createMixingList(myStatus)
77 <    if (myStatus /= 0) then
78 <       status = -1
79 <       return
75 >    havePolicy = .true.
76 >
77 >    if (haveCut) then
78 >       status = 0
79 >       call createMixingList(myStatus)
80 >       if (myStatus /= 0) then
81 >          status = -1
82 >          return
83 >       end if
84 >      
85 >       LJ_FF_initialized = .true.
86      end if
87 <    
81 <    LJ_FF_initialized = .true.
82 <    
87 >  
88    end subroutine init_LJ_FF
89    
90 <  subroutine LJ_new_rcut(rcut, status)
90 >  subroutine setCutoffLJ(rcut, status)
91      integer :: status, myStatus
92      real(kind=dp) :: rcut
93  
94 <    LJ_rcut = rcut
94 > #define __FORTRAN90
95 > #include "fSwitchingFunction.h"
96 >
97      status = 0
91    call createMixingList(myStatus)
92    if (myStatus /= 0) then
93       status = -1
94       return
95    end if
98  
99 +    LJ_rcut = rcut
100 +    call set_switch(LJ_SWITCH, rcut, rcut)
101 +
102 + !!$    ! ATTENTION! This is a hardwiring of rcut!
103 + !!$    LJ_rcut = 9.0d0
104 +    haveCut = .true.
105 +
106 +    if (havePolicy) then
107 +       status = 0
108 +       call createMixingList(myStatus)
109 +       if (myStatus /= 0) then
110 +          status = -1
111 +          return
112 +       end if
113 +      
114 +       LJ_FF_initialized = .true.
115 +    end if    
116 +    
117      return
118 <  end subroutine LJ_new_rcut
118 >  end subroutine setCutoffLJ
119    
120    subroutine createMixingList(status)
121      integer :: nAtypes
# Line 115 | Line 135 | contains
135          
136      if (.not. associated(ljMixed)) then
137         allocate(ljMixed(nAtypes, nAtypes))
138 <    else
119 <       status = -1
120 <       return
121 <    end if
138 >    endif
139  
140      rcut6 = LJ_rcut**6
141  
142 + ! This loops through all atypes, even those that don't support LJ forces.
143      do i = 1, nAtypes
144  
145         call getElementProperty(atypes, i, "lj_epsilon", myEpsilon_i)
# Line 129 | Line 147 | contains
147         ! do self mixing rule
148         ljMixed(i,i)%sigma   = mySigma_i
149        
150 <       ljMixed(i,i)%sigma6  = (ljMixed(i,i)%sigma) ** 6                                                                                                  
133 <       ljMixed(i,i)%tp6     = ljMixed(i,i)%sigma6/rcut6
150 >       ljMixed(i,i)%sigma6  = (ljMixed(i,i)%sigma) ** 6
151  
152 +       ljMixed(i,i)%tp6     = (ljMixed(i,i)%sigma6)/rcut6
153 +
154         ljMixed(i,i)%tp12    = (ljMixed(i,i)%tp6) ** 2
155  
156 +
157         ljMixed(i,i)%epsilon = myEpsilon_i
158  
159         ljMixed(i,i)%delta = -4.0_DP * ljMixed(i,i)%epsilon * &
# Line 176 | Line 196 | contains
196      
197    end subroutine createMixingList
198    
199 <  subroutine do_lj_pair(atom1, atom2, d, rij, r2, pot, f, do_pot, do_stress)
199 >  subroutine do_lj_pair(atom1, atom2, d, rij, r2, sw, vpair, pot, f, &
200 >       do_pot, do_stress)
201  
202      integer, intent(in) ::  atom1, atom2
203      real( kind = dp ), intent(in) :: rij, r2
204 <    real( kind = dp ) :: pot
205 <    real( kind = dp ), dimension(3,getNlocal()) :: f    
204 >    real( kind = dp ) :: pot, sw, vpair
205 >    real( kind = dp ), dimension(3,nLocal) :: f    
206      real( kind = dp ), intent(in), dimension(3) :: d
207      logical, intent(in) :: do_pot, do_stress
208  
# Line 195 | Line 216 | contains
216      real( kind = dp ) :: t6
217      real( kind = dp ) :: t12
218      real( kind = dp ) :: delta
219 +    integer :: id1, id2
220  
199
221      if (rij.lt.LJ_rcut)  then
222  
223         ! Look up the correct parameters in the mixing matrix
224   #ifdef IS_MPI
225         sigma6   = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma6
226         epsilon  = ljMixed(atid_Row(atom1),atid_Col(atom2))%epsilon
227 <       delta    = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta
227 >       !delta    = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta
228   #else
229         sigma6   = ljMixed(atid(atom1),atid(atom2))%sigma6
230         epsilon  = ljMixed(atid(atom1),atid(atom2))%epsilon
231 <       delta    = ljMixed(atid(atom1),atid(atom2))%delta
231 >       !delta    = ljMixed(atid(atom1),atid(atom2))%delta
232   #endif
233        
234         r6 = r2 * r2 * r2
# Line 215 | Line 236 | contains
236         t6  = sigma6/ r6
237         t12 = t6 * t6    
238              
239 <       pot_temp = 4.0E0_DP * epsilon * (t12 - t6) + delta      
239 >       !pot_temp = 4.0E0_DP * epsilon * (t12 - t6) + delta      
240 >       pot_temp = 4.0E0_DP * epsilon * (t12 - t6) * sw
241 >       vpair = vpair + pot_temp
242        
243 <       dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
243 >       dudr = sw * 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
244        
245         drdx = d(1) / rij
246         drdy = d(2) / rij
# Line 227 | Line 250 | contains
250         fy = dudr * drdy
251         fz = dudr * drdz
252        
253 +      
254   #ifdef IS_MPI
255         if (do_pot) then
256            pot_Row(atom1) = pot_Row(atom1) + pot_temp*0.5
# Line 255 | Line 279 | contains
279        
280         if (do_stress) then
281  
282 <          if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then
283 <             tau_Temp(1) = tau_Temp(1) + fx * d(1)
284 <             tau_Temp(2) = tau_Temp(2) + fx * d(2)
285 <             tau_Temp(3) = tau_Temp(3) + fx * d(3)
286 <             tau_Temp(4) = tau_Temp(4) + fy * d(1)
287 <             tau_Temp(5) = tau_Temp(5) + fy * d(2)
288 <             tau_Temp(6) = tau_Temp(6) + fy * d(3)
289 <             tau_Temp(7) = tau_Temp(7) + fz * d(1)
290 <             tau_Temp(8) = tau_Temp(8) + fz * d(2)
291 <             tau_Temp(9) = tau_Temp(9) + fz * d(3)
282 > #ifdef IS_MPI
283 >          id1 = tagRow(atom1)
284 >          id2 = tagColumn(atom2)
285 > #else
286 >          id1 = atom1
287 >          id2 = atom2
288 > #endif
289 >
290 >          if (molMembershipList(id1) .ne. molMembershipList(id2)) then
291 >
292 >             ! because the d vector is the rj - ri vector, and
293 >             ! because fx, fy, fz are the force on atom i, we need a
294 >             ! negative sign here:
295 >
296 >             tau_Temp(1) = tau_Temp(1) - d(1) * fx
297 >             tau_Temp(2) = tau_Temp(2) - d(1) * fy
298 >             tau_Temp(3) = tau_Temp(3) - d(1) * fz
299 >             tau_Temp(4) = tau_Temp(4) - d(2) * fx
300 >             tau_Temp(5) = tau_Temp(5) - d(2) * fy
301 >             tau_Temp(6) = tau_Temp(6) - d(2) * fz
302 >             tau_Temp(7) = tau_Temp(7) - d(3) * fx
303 >             tau_Temp(8) = tau_Temp(8) - d(3) * fy
304 >             tau_Temp(9) = tau_Temp(9) - d(3) * fz
305 >            
306               virial_Temp = virial_Temp + &
307                    (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
270          endif
308  
309 +          endif
310         endif
311        
312      endif
# Line 291 | Line 329 | contains
329      
330      if (present(status)) status = 0
331      select case (LJ_Mixing_Policy)
332 <    case (LB_MIXING_RULE)
332 >    case (1)
333         select case (thisParam)
334         case ("sigma")
335            myMixParam = 0.5_dp * (param1 + param2)

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines