ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_LJ_FF.F90
Revision: 619
Committed: Tue Jul 15 22:22:41 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 8921 byte(s)
Log Message:
fixed a long lived bug in do_forces. Rrf was not being used in the neighborlist correctly. rcut was conssistently being set lowere than Rrf causing the dipole cutoff region to be to small. Also led to the removal of the taper region to buffer the dipole cutoff.

File Contents

# User Rev Content
1 mmeineke 377 !! Calculates Long Range forces Lennard-Jones interactions.
2     !! Corresponds to the force field defined in lj_FF.cpp
3     !! @author Charles F. Vardeman II
4     !! @author Matthew Meineke
5 mmeineke 619 !! @version $Id: calc_LJ_FF.F90,v 1.8 2003-07-15 22:22:41 mmeineke Exp $, $Date: 2003-07-15 22:22:41 $, $Name: not supported by cvs2svn $, $Revision: 1.8 $
6 mmeineke 377
7     module lj
8     use definitions
9     use atype_module
10     use vector_class
11 chuckv 460 use simulation
12 mmeineke 377 #ifdef IS_MPI
13     use mpiSimulation
14     #endif
15     use force_globals
16    
17     implicit none
18     PRIVATE
19    
20     #define __FORTRAN90
21     #include "fForceField.h"
22    
23     integer, save :: LJ_Mixing_Policy
24 mmeineke 619 real(kind=DP), save :: LJ_rcut
25 mmeineke 377
26     !! Logical has lj force field module been initialized?
27    
28     logical, save :: LJ_FF_initialized = .false.
29    
30     !! Public methods and data
31     public :: init_LJ_FF
32     public :: LJ_new_rcut
33     public :: do_lj_pair
34    
35     type :: lj_mixed_params
36     !! Lennard-Jones epsilon
37     real ( kind = dp ) :: epsilon = 0.0_dp
38     !! Lennard-Jones Sigma
39     real ( kind = dp ) :: sigma = 0.0_dp
40     !! Lennard-Jones Sigma to sixth
41     real ( kind = dp ) :: sigma6 = 0.0_dp
42     !!
43     real ( kind = dp ) :: tp6
44     real ( kind = dp ) :: tp12
45    
46     real ( kind = dp ) :: delta = 0.0_dp
47    
48     end type lj_mixed_params
49    
50     type (lj_mixed_params), dimension(:,:), pointer :: ljMixed
51    
52     contains
53    
54     subroutine init_LJ_FF(mix_Policy, rcut, status)
55     integer, intent(in) :: mix_Policy
56     integer, intent(out) :: status
57     integer :: myStatus
58     real(kind=dp) :: rcut
59    
60     if (mix_Policy == LB_MIXING_RULE) then
61     LJ_Mixing_Policy = LB_MIXING_RULE
62     else
63     if (mix_Policy == EXPLICIT_MIXING_RULE) then
64     LJ_Mixing_Policy = EXPLICIT_MIXING_RULE
65     else
66     write(*,*) 'Unknown Mixing Policy!'
67     status = -1
68     return
69     endif
70     endif
71    
72     LJ_rcut = rcut
73    
74     status = 0
75     call createMixingList(myStatus)
76     if (myStatus /= 0) then
77     status = -1
78     return
79     end if
80    
81     LJ_FF_initialized = .true.
82    
83     end subroutine init_LJ_FF
84    
85     subroutine LJ_new_rcut(rcut, status)
86     integer :: status, myStatus
87     real(kind=dp) :: rcut
88 mmeineke 619
89 mmeineke 377 status = 0
90    
91 mmeineke 619 if ( rcut .ne. LJ_rcut ) then
92     LJ_rcut = rcut
93     call createMixingList(myStatus)
94     if (myStatus /= 0) then
95     status = -1
96     return
97     end if
98     endif
99    
100    
101 mmeineke 377 return
102     end subroutine LJ_new_rcut
103    
104     subroutine createMixingList(status)
105     integer :: nAtypes
106     integer :: status
107     integer :: i
108     integer :: j
109     real ( kind = dp ) :: mySigma_i,mySigma_j
110     real ( kind = dp ) :: myEpsilon_i,myEpsilon_j
111     real ( kind = dp ) :: rcut6
112     status = 0
113    
114     nAtypes = getSize(atypes)
115     if (nAtypes == 0) then
116     status = -1
117     return
118     end if
119    
120     if (.not. associated(ljMixed)) then
121     allocate(ljMixed(nAtypes, nAtypes))
122     else
123     status = -1
124     return
125     end if
126    
127     rcut6 = LJ_rcut**6
128    
129     do i = 1, nAtypes
130    
131     call getElementProperty(atypes, i, "lj_epsilon", myEpsilon_i)
132     call getElementProperty(atypes, i, "lj_sigma", mySigma_i)
133     ! do self mixing rule
134     ljMixed(i,i)%sigma = mySigma_i
135    
136     ljMixed(i,i)%sigma6 = (ljMixed(i,i)%sigma) ** 6
137     ljMixed(i,i)%tp6 = ljMixed(i,i)%sigma6/rcut6
138    
139     ljMixed(i,i)%tp12 = (ljMixed(i,i)%tp6) ** 2
140    
141     ljMixed(i,i)%epsilon = myEpsilon_i
142    
143     ljMixed(i,i)%delta = -4.0_DP * ljMixed(i,i)%epsilon * &
144     (ljMixed(i,i)%tp12 - ljMixed(i,i)%tp6)
145    
146     do j = i + 1, nAtypes
147     call getElementProperty(atypes,j,"lj_epsilon",myEpsilon_j)
148     call getElementProperty(atypes,j,"lj_sigma", mySigma_j)
149    
150     ljMixed(i,j)%sigma = &
151     calcLJMix("sigma",mySigma_i, &
152     mySigma_j)
153    
154     ljMixed(i,j)%sigma6 = &
155     (ljMixed(i,j)%sigma)**6
156    
157    
158     ljMixed(i,j)%tp6 = ljMixed(i,j)%sigma6/rcut6
159    
160     ljMixed(i,j)%tp12 = (ljMixed(i,j)%tp6) ** 2
161    
162    
163     ljMixed(i,j)%epsilon = &
164     calcLJMix("epsilon",myEpsilon_i, &
165     myEpsilon_j)
166    
167     ljMixed(i,j)%delta = -4.0_DP * ljMixed(i,j)%epsilon * &
168     (ljMixed(i,j)%tp12 - ljMixed(i,j)%tp6)
169    
170    
171     ljMixed(j,i)%sigma = ljMixed(i,j)%sigma
172     ljMixed(j,i)%sigma6 = ljMixed(i,j)%sigma6
173     ljMixed(j,i)%tp6 = ljMixed(i,j)%tp6
174     ljMixed(j,i)%tp12 = ljMixed(i,j)%tp12
175     ljMixed(j,i)%epsilon = ljMixed(i,j)%epsilon
176     ljMixed(j,i)%delta = ljMixed(i,j)%delta
177    
178     end do
179     end do
180    
181     end subroutine createMixingList
182    
183     subroutine do_lj_pair(atom1, atom2, d, rij, r2, pot, f, do_pot, do_stress)
184    
185     integer, intent(in) :: atom1, atom2
186     real( kind = dp ), intent(in) :: rij, r2
187     real( kind = dp ) :: pot
188 chuckv 460 real( kind = dp ), dimension(3,getNlocal()) :: f
189 mmeineke 377 real( kind = dp ), intent(in), dimension(3) :: d
190     logical, intent(in) :: do_pot, do_stress
191    
192     ! local Variables
193     real( kind = dp ) :: drdx, drdy, drdz
194     real( kind = dp ) :: fx, fy, fz
195     real( kind = dp ) :: pot_temp, dudr
196     real( kind = dp ) :: sigma6
197     real( kind = dp ) :: epsilon
198     real( kind = dp ) :: r6
199     real( kind = dp ) :: t6
200     real( kind = dp ) :: t12
201     real( kind = dp ) :: delta
202 gezelter 490 integer :: id1, id2
203 mmeineke 377
204    
205     if (rij.lt.LJ_rcut) then
206    
207     ! Look up the correct parameters in the mixing matrix
208     #ifdef IS_MPI
209     sigma6 = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma6
210     epsilon = ljMixed(atid_Row(atom1),atid_Col(atom2))%epsilon
211     delta = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta
212     #else
213     sigma6 = ljMixed(atid(atom1),atid(atom2))%sigma6
214     epsilon = ljMixed(atid(atom1),atid(atom2))%epsilon
215     delta = ljMixed(atid(atom1),atid(atom2))%delta
216     #endif
217    
218     r6 = r2 * r2 * r2
219    
220     t6 = sigma6/ r6
221     t12 = t6 * t6
222    
223     pot_temp = 4.0E0_DP * epsilon * (t12 - t6) + delta
224    
225     dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
226    
227 chuckv 482 drdx = d(1) / rij
228     drdy = d(2) / rij
229     drdz = d(3) / rij
230 mmeineke 377
231     fx = dudr * drdx
232     fy = dudr * drdy
233     fz = dudr * drdz
234    
235     #ifdef IS_MPI
236     if (do_pot) then
237     pot_Row(atom1) = pot_Row(atom1) + pot_temp*0.5
238     pot_Col(atom2) = pot_Col(atom2) + pot_temp*0.5
239     endif
240    
241     f_Row(1,atom1) = f_Row(1,atom1) + fx
242     f_Row(2,atom1) = f_Row(2,atom1) + fy
243     f_Row(3,atom1) = f_Row(3,atom1) + fz
244    
245     f_Col(1,atom2) = f_Col(1,atom2) - fx
246     f_Col(2,atom2) = f_Col(2,atom2) - fy
247     f_Col(3,atom2) = f_Col(3,atom2) - fz
248    
249     #else
250     if (do_pot) pot = pot + pot_temp
251    
252     f(1,atom1) = f(1,atom1) + fx
253     f(2,atom1) = f(2,atom1) + fy
254     f(3,atom1) = f(3,atom1) + fz
255    
256     f(1,atom2) = f(1,atom2) - fx
257     f(2,atom2) = f(2,atom2) - fy
258     f(3,atom2) = f(3,atom2) - fz
259     #endif
260    
261     if (do_stress) then
262 gezelter 483
263 gezelter 490 #ifdef IS_MPI
264     id1 = tagRow(atom1)
265     id2 = tagColumn(atom2)
266     #else
267     id1 = atom1
268     id2 = atom2
269     #endif
270    
271     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
272 gezelter 611
273     ! because the d vector is the rj - ri vector, and
274     ! because fx, fy, fz are the force on atom i, we need a
275     ! negative sign here:
276 mmeineke 491
277 gezelter 611 tau_Temp(1) = tau_Temp(1) - d(1) * fx
278     tau_Temp(2) = tau_Temp(2) - d(1) * fy
279     tau_Temp(3) = tau_Temp(3) - d(1) * fz
280     tau_Temp(4) = tau_Temp(4) - d(2) * fx
281     tau_Temp(5) = tau_Temp(5) - d(2) * fy
282     tau_Temp(6) = tau_Temp(6) - d(2) * fz
283     tau_Temp(7) = tau_Temp(7) - d(3) * fx
284     tau_Temp(8) = tau_Temp(8) - d(3) * fy
285     tau_Temp(9) = tau_Temp(9) - d(3) * fz
286    
287 gezelter 483 virial_Temp = virial_Temp + &
288     (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
289 mmeineke 491
290 gezelter 483 endif
291 mmeineke 377 endif
292    
293     endif
294     return
295    
296     end subroutine do_lj_pair
297    
298    
299     !! Calculates the mixing for sigma or epslon
300    
301     function calcLJMix(thisParam,param1,param2,status) result(myMixParam)
302     character(len=*) :: thisParam
303     real(kind = dp) :: param1
304     real(kind = dp) :: param2
305     real(kind = dp ) :: myMixParam
306    
307     integer, optional :: status
308    
309     myMixParam = 0.0_dp
310    
311     if (present(status)) status = 0
312     select case (LJ_Mixing_Policy)
313     case (LB_MIXING_RULE)
314     select case (thisParam)
315     case ("sigma")
316     myMixParam = 0.5_dp * (param1 + param2)
317     case ("epsilon")
318     myMixParam = sqrt(param1 * param2)
319     case default
320     status = -1
321     end select
322     case default
323     status = -1
324     end select
325     end function calcLJMix
326    
327     end module lj