ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_LJ_FF.F90
Revision: 626
Committed: Wed Jul 16 21:30:56 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 9035 byte(s)
Log Message:
Changed how cutoffs were handled from C. Now notifyCutoffs in Fortran notifies those who need the information of any changes to cutoffs.

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