ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_LJ_FF.F90
Revision: 1198
Committed: Thu May 27 00:48:12 2004 UTC (20 years, 1 month ago) by tim
File size: 8551 byte(s)
Log Message:
in the progress of fixing MPI version of cutoff group

File Contents

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