ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_LJ_FF.F90
Revision: 900
Committed: Tue Jan 6 18:54:57 2004 UTC (20 years, 6 months ago) by chuckv
File size: 8995 byte(s)
Log Message:
Making do_Forces a little more sane

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