ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/calc_LJ_FF.F90
Revision: 306
Committed: Mon Mar 10 19:26:45 2003 UTC (21 years, 3 months ago) by chuckv
File size: 7153 byte(s)
Log Message:
More changes to code. Hopefully these will commit smoothly.

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.2 2003-03-10 19:26:45 chuckv Exp $, $Date: 2003-03-10 19:26:45 $, $Name: not supported by cvs2svn $, $Revision: 1.2 $
6
7
8
9 module lj_ff
10 use simulation
11 use definitions
12 use generic_lists
13
14 #ifdef IS_MPI
15 use mpiSimulation
16 #endif
17 implicit none
18 PRIVATE
19
20 !! Logical has lj force field module been initialized?
21 logical, save :: isljFFinit = .false.
22
23
24 !! Public methods and data
25 public :: getLjPot
26 public :: init_ljFF
27
28 type :: lj_mixed_params
29 !! Mass of Particle
30 real ( kind = dp ) :: mass = 0.0_dp
31 !! Lennard-Jones epslon
32 real ( kind = dp ) :: epsilon = 0.0_dp
33 !! Lennard-Jones Sigma
34 real ( kind = dp ) :: sigma = 0.0_dp
35 !! Lennard-Jones Sigma Squared
36 real ( kind = dp ) :: sigma2 = 0.0_dp
37 !! Lennard-Jones Sigma to sixth
38 real ( kind = dp ) :: sigma6 = 0.0_dp
39 end type lj_mixed_params
40
41 type (lj_mixed_params), dimension(:,:) :: ljMixed
42
43 contains
44
45 subroutine init_ljFF(ListHead,status)
46 type (atype),pointer :: ListHead
47 integer :: nTypes
48 integer, intent(out) :: status
49
50 integer :: myStatus
51
52 status = 0
53 call createMixingList(ListHead,myStatus)
54 if (myStatus /= 0) then
55 status = -1
56 return
57 end if
58
59 end subroutine init_ljFF
60
61 subroutine createMixingList(ListHead,status)
62 type (atype), pointer :: ListHead
63 integer :: listSize
64 integer :: status
65 integer :: i
66 integer :: j
67
68 integer :: i = 0
69 integer :: j = 0
70 type (atype), pointer :: tmpPtr_i => null()
71 type (atype), pointer :: tmpPtr_j => null()
72 status = 0
73
74 listSize = getListLen(ListHead)
75 if (listSize == 0) then
76 status = -1
77 return
78 end if
79
80
81 if (.not. associated(ljMixed)) then
82 allocate(ljMixed(listSize,listSize))
83 else
84 status = -1
85 return
86 end if
87
88
89
90 tmpPtr_i => ListHead
91 tmpPtr_j => tmpPtr_i%next
92 do while (associated(tmpPtr_i))
93 i = i + 1
94 ! do self mixing rule
95 ljMixed(i,i)%sigma = tmpPtr_i%sigma
96
97 ljMixed(i,i)%sigma2 = (ljMixed(i,i)%sigma) ** 2
98
99 ljMixed(i,i)%sigma6 = (ljMixed(i,i)%sigma) ** 6
100 ljMixed(i,i)%epsilon = tmpPtr_i%epsilon
101
102 j = i + 1
103 do while (associated(tmpPtr_j))
104
105 ljMixed(i,j)%sigma = &
106 calcLJMix("sigma",tmpPtr_i%sigma, &
107 tmpPtr_j%sigma)
108
109 ljMixed(i,j)%sigma2 = &
110 (ljMixed(i,j)%sigma)**2
111
112 ljMixed(i,j)%sigma6 = &
113 (ljMixed(i,j)%sigma)**6
114
115 ljMixed(i,j)%epsilon = &
116 calcLJMix("epsilon",tmpPtr_i%epsilon, &
117 tmpPtr_j%epsilon)
118 ljMixed(j,i)%sigma = ljMixed(i,j)%sigma
119 ljMixed(j,i)%sigma2 = ljMixed(i,j)%sigma2
120 ljMixed(j,i)%sigma6 = ljMixed(i,j)%sigma6
121 ljMixed(j,i)%epsilon = ljMixed(i,j)%epsilon
122
123
124 tmpPtr_j => tmpPtr_j%next
125 j = j + 1
126 end do
127 ! advance pointers
128 tmpPtr_i => tmpPtr_i%next
129 if (associated(tmpPtr_i)) then
130 tmpPtr_j => tmpPtr_i%next
131 endif
132
133 end do
134
135 end subroutine createMixingList
136
137
138
139
140
141 !! Calculates the potential between two lj particles based on two lj_atype pointers, optionally returns second
142 !! derivatives.
143 subroutine getLjForce(r,pot,dudr,atype1,atype2,d2,status)
144 ! arguments
145 !! Length of vector between particles
146 real( kind = dp ), intent(in) :: r
147 !! Potential Energy
148 real( kind = dp ), intent(out) :: pot
149 !! Derivatve wrt postion
150 real( kind = dp ), intent(out) :: dudr
151 !! Second Derivative, optional, used mainly for normal mode calculations.
152 real( kind = dp ), intent(out), optional :: d2
153
154 type (atype), pointer :: atype1
155 type (atype), pointer :: atype2
156
157 integer, intent(out), optional :: status
158
159 ! local Variables
160 real( kind = dp ) :: sigma
161 real( kind = dp ) :: sigma2
162 real( kind = dp ) :: sigma6
163 real( kind = dp ) :: epsilon
164
165 real( kind = dp ) :: rcut
166 real( kind = dp ) :: rcut2
167 real( kind = dp ) :: rcut6
168 real( kind = dp ) :: r2
169 real( kind = dp ) :: r6
170
171 real( kind = dp ) :: t6
172 real( kind = dp ) :: t12
173 real( kind = dp ) :: tp6
174 real( kind = dp ) :: tp12
175 real( kind = dp ) :: delta
176
177 logical :: doSec = .false.
178
179 integer :: errorStat
180
181 !! Optional Argument Checking
182 ! Check to see if we need to do second derivatives
183
184 if (present(d2)) doSec = .true.
185 if (present(status)) status = 0
186
187 ! Look up the correct parameters in the mixing matrix
188 sigma = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma
189 sigma2 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma2
190 sigma6 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma6
191 epsilon = ljMixed(atype1%atype_ident,atype2%atype_ident)%epsilon
192
193
194
195
196 call getRcut(rcut,rcut2=rcut2,rcut6=rcut6,status=errorStat)
197
198 r2 = r * r
199 r6 = r2 * r2 * r2
200
201 t6 = sigma6/ r6
202 t12 = t6 * t6
203
204
205
206 tp6 = sigma6 / rcut6
207 tp12 = tp6*tp6
208
209 delta = -4.0E0_DP*epsilon * (tp12 - tp6)
210
211 if (r.le.rcut) then
212 pot = 4.0E0_DP * epsilon * (t12 - t6) + delta
213 dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / r
214 if(doSec) d2 = 24.0E0_DP * epsilon * (26.0E0_DP*t12 - 7.0E0_DP*t6)/r/r
215 else
216 pot = 0.0E0_DP
217 dudr = 0.0E0_DP
218 if(doSec) d2 = 0.0E0_DP
219 endif
220
221 return
222
223
224
225 end subroutine getLJForce
226
227
228 !! Calculates the mixing for sigma or epslon based on character string for initialzition of mixing array.
229 function calcLJMix(thisParam,param1,param2,status) result(myMixParam)
230 character(len=*) :: thisParam
231 real(kind = dp) :: param1
232 real(kind = dp) :: param2
233 real(kind = dp ) :: myMixParam
234 character(len = getStringLen()) :: thisMixingRule
235 integer, optional :: status
236
237 !! get the mixing rules from the simulation
238 thisMixingRule = returnMixingRules()
239 myMixParam = 0.0_dp
240
241 if (present(status)) status = 0
242 select case (thisMixingRule)
243 case ("standard")
244 select case (thisParam)
245 case ("sigma")
246 myMixParam = 0.5_dp * (param1 + param2)
247 case ("epsilon")
248 myMixParam = sqrt(param1 * param2)
249 case default
250 status = -1
251 end select
252 case("LJglass")
253 case default
254 status = -1
255 end select
256 end function calcLJMix
257
258
259
260 end module lj_ff