ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_LJ_FF.F90
Revision: 623
Committed: Wed Jul 16 16:40:03 2003 UTC (20 years, 11 months ago) by chuckv
File size: 8876 byte(s)
Log Message:
Fixed bug in updating mixing lists

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 chuckv 623 !! @version $Id: calc_LJ_FF.F90,v 1.9 2003-07-16 16:40:03 chuckv Exp $, $Date: 2003-07-16 16:40:03 $, $Name: not supported by cvs2svn $, $Revision: 1.9 $
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 chuckv 623 endif
123 mmeineke 377
124     rcut6 = LJ_rcut**6
125    
126     do i = 1, nAtypes
127    
128     call getElementProperty(atypes, i, "lj_epsilon", myEpsilon_i)
129     call getElementProperty(atypes, i, "lj_sigma", mySigma_i)
130     ! do self mixing rule
131     ljMixed(i,i)%sigma = mySigma_i
132    
133     ljMixed(i,i)%sigma6 = (ljMixed(i,i)%sigma) ** 6
134     ljMixed(i,i)%tp6 = ljMixed(i,i)%sigma6/rcut6
135    
136     ljMixed(i,i)%tp12 = (ljMixed(i,i)%tp6) ** 2
137    
138     ljMixed(i,i)%epsilon = myEpsilon_i
139    
140     ljMixed(i,i)%delta = -4.0_DP * ljMixed(i,i)%epsilon * &
141     (ljMixed(i,i)%tp12 - ljMixed(i,i)%tp6)
142    
143     do j = i + 1, nAtypes
144     call getElementProperty(atypes,j,"lj_epsilon",myEpsilon_j)
145     call getElementProperty(atypes,j,"lj_sigma", mySigma_j)
146    
147     ljMixed(i,j)%sigma = &
148     calcLJMix("sigma",mySigma_i, &
149     mySigma_j)
150    
151     ljMixed(i,j)%sigma6 = &
152     (ljMixed(i,j)%sigma)**6
153    
154    
155     ljMixed(i,j)%tp6 = ljMixed(i,j)%sigma6/rcut6
156    
157     ljMixed(i,j)%tp12 = (ljMixed(i,j)%tp6) ** 2
158    
159    
160     ljMixed(i,j)%epsilon = &
161     calcLJMix("epsilon",myEpsilon_i, &
162     myEpsilon_j)
163    
164     ljMixed(i,j)%delta = -4.0_DP * ljMixed(i,j)%epsilon * &
165     (ljMixed(i,j)%tp12 - ljMixed(i,j)%tp6)
166    
167    
168     ljMixed(j,i)%sigma = ljMixed(i,j)%sigma
169     ljMixed(j,i)%sigma6 = ljMixed(i,j)%sigma6
170     ljMixed(j,i)%tp6 = ljMixed(i,j)%tp6
171     ljMixed(j,i)%tp12 = ljMixed(i,j)%tp12
172     ljMixed(j,i)%epsilon = ljMixed(i,j)%epsilon
173     ljMixed(j,i)%delta = ljMixed(i,j)%delta
174    
175     end do
176     end do
177    
178     end subroutine createMixingList
179    
180     subroutine do_lj_pair(atom1, atom2, d, rij, r2, pot, f, do_pot, do_stress)
181    
182     integer, intent(in) :: atom1, atom2
183     real( kind = dp ), intent(in) :: rij, r2
184     real( kind = dp ) :: pot
185 chuckv 460 real( kind = dp ), dimension(3,getNlocal()) :: f
186 mmeineke 377 real( kind = dp ), intent(in), dimension(3) :: d
187     logical, intent(in) :: do_pot, do_stress
188    
189     ! local Variables
190     real( kind = dp ) :: drdx, drdy, drdz
191     real( kind = dp ) :: fx, fy, fz
192     real( kind = dp ) :: pot_temp, dudr
193     real( kind = dp ) :: sigma6
194     real( kind = dp ) :: epsilon
195     real( kind = dp ) :: r6
196     real( kind = dp ) :: t6
197     real( kind = dp ) :: t12
198     real( kind = dp ) :: delta
199 gezelter 490 integer :: id1, id2
200 mmeineke 377
201    
202     if (rij.lt.LJ_rcut) then
203    
204     ! Look up the correct parameters in the mixing matrix
205     #ifdef IS_MPI
206     sigma6 = ljMixed(atid_Row(atom1),atid_Col(atom2))%sigma6
207     epsilon = ljMixed(atid_Row(atom1),atid_Col(atom2))%epsilon
208     delta = ljMixed(atid_Row(atom1),atid_Col(atom2))%delta
209     #else
210     sigma6 = ljMixed(atid(atom1),atid(atom2))%sigma6
211     epsilon = ljMixed(atid(atom1),atid(atom2))%epsilon
212     delta = ljMixed(atid(atom1),atid(atom2))%delta
213     #endif
214    
215     r6 = r2 * r2 * r2
216    
217     t6 = sigma6/ r6
218     t12 = t6 * t6
219    
220     pot_temp = 4.0E0_DP * epsilon * (t12 - t6) + delta
221    
222     dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / rij
223    
224 chuckv 482 drdx = d(1) / rij
225     drdy = d(2) / rij
226     drdz = d(3) / rij
227 mmeineke 377
228     fx = dudr * drdx
229     fy = dudr * drdy
230     fz = dudr * drdz
231    
232     #ifdef IS_MPI
233     if (do_pot) then
234     pot_Row(atom1) = pot_Row(atom1) + pot_temp*0.5
235     pot_Col(atom2) = pot_Col(atom2) + pot_temp*0.5
236     endif
237    
238     f_Row(1,atom1) = f_Row(1,atom1) + fx
239     f_Row(2,atom1) = f_Row(2,atom1) + fy
240     f_Row(3,atom1) = f_Row(3,atom1) + fz
241    
242     f_Col(1,atom2) = f_Col(1,atom2) - fx
243     f_Col(2,atom2) = f_Col(2,atom2) - fy
244     f_Col(3,atom2) = f_Col(3,atom2) - fz
245    
246     #else
247     if (do_pot) pot = pot + pot_temp
248    
249     f(1,atom1) = f(1,atom1) + fx
250     f(2,atom1) = f(2,atom1) + fy
251     f(3,atom1) = f(3,atom1) + fz
252    
253     f(1,atom2) = f(1,atom2) - fx
254     f(2,atom2) = f(2,atom2) - fy
255     f(3,atom2) = f(3,atom2) - fz
256     #endif
257    
258     if (do_stress) then
259 gezelter 483
260 gezelter 490 #ifdef IS_MPI
261     id1 = tagRow(atom1)
262     id2 = tagColumn(atom2)
263     #else
264     id1 = atom1
265     id2 = atom2
266     #endif
267    
268     if (molMembershipList(id1) .ne. molMembershipList(id2)) then
269 gezelter 611
270     ! because the d vector is the rj - ri vector, and
271     ! because fx, fy, fz are the force on atom i, we need a
272     ! negative sign here:
273 mmeineke 491
274 gezelter 611 tau_Temp(1) = tau_Temp(1) - d(1) * fx
275     tau_Temp(2) = tau_Temp(2) - d(1) * fy
276     tau_Temp(3) = tau_Temp(3) - d(1) * fz
277     tau_Temp(4) = tau_Temp(4) - d(2) * fx
278     tau_Temp(5) = tau_Temp(5) - d(2) * fy
279     tau_Temp(6) = tau_Temp(6) - d(2) * fz
280     tau_Temp(7) = tau_Temp(7) - d(3) * fx
281     tau_Temp(8) = tau_Temp(8) - d(3) * fy
282     tau_Temp(9) = tau_Temp(9) - d(3) * fz
283    
284 gezelter 483 virial_Temp = virial_Temp + &
285     (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
286 mmeineke 491
287 gezelter 483 endif
288 mmeineke 377 endif
289    
290     endif
291     return
292    
293     end subroutine do_lj_pair
294    
295    
296     !! Calculates the mixing for sigma or epslon
297    
298     function calcLJMix(thisParam,param1,param2,status) result(myMixParam)
299     character(len=*) :: thisParam
300     real(kind = dp) :: param1
301     real(kind = dp) :: param2
302     real(kind = dp ) :: myMixParam
303    
304     integer, optional :: status
305    
306     myMixParam = 0.0_dp
307    
308     if (present(status)) status = 0
309     select case (LJ_Mixing_Policy)
310     case (LB_MIXING_RULE)
311     select case (thisParam)
312     case ("sigma")
313     myMixParam = 0.5_dp * (param1 + param2)
314     case ("epsilon")
315     myMixParam = sqrt(param1 * param2)
316     case default
317     status = -1
318     end select
319     case default
320     status = -1
321     end select
322     end function calcLJMix
323    
324     end module lj