1 |
chuckv |
239 |
!! 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 |
254 |
!! @version $Id: lj_FF.F90,v 1.14 2003-01-30 20:03:36 chuckv Exp $, $Date: 2003-01-30 20:03:36 $, $Name: not supported by cvs2svn $, $Revision: 1.14 $ |
6 |
chuckv |
239 |
|
7 |
|
|
|
8 |
|
|
|
9 |
chuckv |
215 |
module lj_ff |
10 |
|
|
use simulation |
11 |
chuckv |
246 |
use definitions |
12 |
chuckv |
247 |
use generic_atypes |
13 |
chuckv |
230 |
#ifdef IS_MPI |
14 |
|
|
use mpiSimulation |
15 |
|
|
#endif |
16 |
chuckv |
215 |
implicit none |
17 |
chuckv |
230 |
PRIVATE |
18 |
chuckv |
215 |
|
19 |
chuckv |
230 |
!! Number of lj_atypes in lj_atype_list |
20 |
chuckv |
215 |
integer, save :: n_lj_atypes = 0 |
21 |
|
|
|
22 |
chuckv |
247 |
!! Global list of lj atypes in simulation |
23 |
|
|
type (lj_atype), pointer :: ljListHead => null() |
24 |
|
|
type (lj_atype), pointer :: ljListTail => null() |
25 |
chuckv |
230 |
|
26 |
chuckv |
215 |
|
27 |
chuckv |
230 |
!! LJ mixing array |
28 |
chuckv |
246 |
type (lj_atype), dimension(:,:), pointer :: ljMixed => null() |
29 |
chuckv |
215 |
|
30 |
chuckv |
230 |
|
31 |
|
|
!! Neighbor list and commom arrays |
32 |
chuckv |
239 |
!! This should eventally get moved to a neighbor list type |
33 |
chuckv |
230 |
integer, allocatable, dimension(:) :: point |
34 |
|
|
integer, allocatable, dimension(:) :: list |
35 |
chuckv |
239 |
integer, parameter :: listMultiplier = 80 |
36 |
|
|
integer :: nListAllocs = 0 |
37 |
|
|
integer, parameter :: maxListAllocs = 5 |
38 |
chuckv |
230 |
|
39 |
chuckv |
253 |
logical, save :: firstTime = .True. |
40 |
chuckv |
247 |
|
41 |
|
|
!! Atype identity pointer lists |
42 |
chuckv |
230 |
#ifdef IS_MPI |
43 |
chuckv |
246 |
!! Row lj_atype pointer list |
44 |
chuckv |
253 |
type (lj_identPtrList), dimension(:), pointer :: identPtrListRow => null() |
45 |
chuckv |
246 |
!! Column lj_atype pointer list |
46 |
chuckv |
253 |
type (lj_identPtrList), dimension(:), pointer :: identPtrListColumn => null() |
47 |
chuckv |
247 |
#else |
48 |
|
|
type( lj_identPtrList ), dimension(:), pointer :: identPtrList => null() |
49 |
chuckv |
230 |
#endif |
50 |
|
|
|
51 |
|
|
|
52 |
chuckv |
246 |
!! Logical has lj force field module been initialized? |
53 |
|
|
logical, save :: isljFFinit = .false. |
54 |
chuckv |
230 |
|
55 |
|
|
|
56 |
|
|
!! Public methods and data |
57 |
chuckv |
215 |
public :: new_lj_atype |
58 |
chuckv |
230 |
public :: do_lj_ff |
59 |
|
|
public :: getLjPot |
60 |
chuckv |
246 |
public :: init_ljFF |
61 |
chuckv |
215 |
|
62 |
chuckv |
230 |
|
63 |
|
|
|
64 |
|
|
|
65 |
chuckv |
215 |
contains |
66 |
|
|
|
67 |
chuckv |
247 |
!! Adds a new lj_atype to the list. |
68 |
chuckv |
246 |
subroutine new_lj_atype(ident,mass,epsilon,sigma,status) |
69 |
chuckv |
215 |
real( kind = dp ), intent(in) :: mass |
70 |
chuckv |
246 |
real( kind = dp ), intent(in) :: epsilon |
71 |
chuckv |
230 |
real( kind = dp ), intent(in) :: sigma |
72 |
|
|
integer, intent(in) :: ident |
73 |
chuckv |
215 |
integer, intent(out) :: status |
74 |
|
|
|
75 |
chuckv |
247 |
type (lj_atype), pointer :: newLJ_atype |
76 |
chuckv |
216 |
integer :: alloc_error |
77 |
|
|
integer :: atype_counter = 0 |
78 |
chuckv |
230 |
integer :: alloc_size |
79 |
chuckv |
247 |
integer :: err_stat |
80 |
chuckv |
215 |
status = 0 |
81 |
|
|
|
82 |
chuckv |
230 |
|
83 |
|
|
|
84 |
|
|
! allocate a new atype |
85 |
chuckv |
247 |
allocate(newLJ_atype,stat=alloc_error) |
86 |
chuckv |
216 |
if (alloc_error /= 0 ) then |
87 |
|
|
status = -1 |
88 |
|
|
return |
89 |
|
|
end if |
90 |
|
|
|
91 |
|
|
! assign our new lj_atype information |
92 |
chuckv |
247 |
newLJ_atype%mass = mass |
93 |
|
|
newLJ_atype%epsilon = epsilon |
94 |
|
|
newLJ_atype%sigma = sigma |
95 |
|
|
newLJ_atype%sigma2 = sigma * sigma |
96 |
|
|
newLJ_atype%sigma6 = newLJ_atype%sigma2 * newLJ_atype%sigma2 & |
97 |
|
|
* newLJ_atype%sigma2 |
98 |
chuckv |
230 |
! assume that this atype will be successfully added |
99 |
chuckv |
247 |
newLJ_atype%atype_ident = ident |
100 |
|
|
newLJ_atype%atype_number = n_lj_atypes + 1 |
101 |
chuckv |
216 |
|
102 |
chuckv |
247 |
call add_atype(newLJ_atype,ljListHead,ljListTail,err_stat) |
103 |
|
|
if (err_stat /= 0 ) then |
104 |
|
|
status = -1 |
105 |
|
|
return |
106 |
chuckv |
230 |
endif |
107 |
|
|
|
108 |
|
|
n_lj_atypes = n_lj_atypes + 1 |
109 |
|
|
|
110 |
|
|
|
111 |
chuckv |
215 |
end subroutine new_lj_atype |
112 |
|
|
|
113 |
|
|
|
114 |
chuckv |
239 |
subroutine init_ljFF(nComponents,ident, status) |
115 |
|
|
!! Number of components in ident array |
116 |
|
|
integer, intent(inout) :: nComponents |
117 |
|
|
!! Array of identities nComponents long corresponding to |
118 |
|
|
!! ljatype ident. |
119 |
|
|
integer, dimension(nComponents),intent(inout) :: ident |
120 |
|
|
!! Result status, success = 0, error = -1 |
121 |
|
|
integer, intent(out) :: Status |
122 |
chuckv |
215 |
|
123 |
chuckv |
246 |
integer :: alloc_stat |
124 |
|
|
|
125 |
chuckv |
239 |
integer :: thisStat |
126 |
chuckv |
246 |
integer :: i |
127 |
chuckv |
230 |
#ifdef IS_MPI |
128 |
chuckv |
239 |
integer, allocatable, dimension(:) :: identRow |
129 |
|
|
integer, allocatable, dimension(:) :: identCol |
130 |
|
|
integer :: nrow |
131 |
|
|
integer :: ncol |
132 |
|
|
#endif |
133 |
|
|
status = 0 |
134 |
chuckv |
254 |
|
135 |
|
|
|
136 |
|
|
|
137 |
|
|
|
138 |
chuckv |
239 |
!! if were're not in MPI, we just update ljatypePtrList |
139 |
|
|
#ifndef IS_MPI |
140 |
chuckv |
247 |
call create_IdentPtrlst(ident,ljListHead,identPtrList,thisStat) |
141 |
chuckv |
239 |
if ( thisStat /= 0 ) then |
142 |
|
|
status = -1 |
143 |
|
|
return |
144 |
|
|
endif |
145 |
chuckv |
247 |
|
146 |
chuckv |
239 |
!! Allocate pointer lists |
147 |
|
|
allocate(point(nComponents),stat=alloc_stat) |
148 |
|
|
if (alloc_stat /=0) then |
149 |
|
|
status = -1 |
150 |
|
|
return |
151 |
|
|
endif |
152 |
|
|
|
153 |
|
|
allocate(list(nComponents*listMultiplier),stat=alloc_stat) |
154 |
|
|
if (alloc_stat /=0) then |
155 |
|
|
status = -1 |
156 |
|
|
return |
157 |
|
|
endif |
158 |
chuckv |
230 |
|
159 |
chuckv |
239 |
! if were're in MPI, we also have to worry about row and col lists |
160 |
|
|
#else |
161 |
chuckv |
254 |
|
162 |
chuckv |
239 |
! We can only set up forces if mpiSimulation has been setup. |
163 |
|
|
if (.not. isMPISimSet()) then |
164 |
chuckv |
254 |
write(default_error,*) "MPI is not set" |
165 |
chuckv |
239 |
status = -1 |
166 |
|
|
return |
167 |
|
|
endif |
168 |
chuckv |
254 |
nrow = getNrow(plan_row) |
169 |
|
|
ncol = getNcol(plan_col) |
170 |
chuckv |
239 |
!! Allocate temperary arrays to hold gather information |
171 |
|
|
allocate(identRow(nrow),stat=alloc_stat) |
172 |
|
|
if (alloc_stat /= 0 ) then |
173 |
|
|
status = -1 |
174 |
|
|
return |
175 |
|
|
endif |
176 |
chuckv |
215 |
|
177 |
chuckv |
239 |
allocate(identCol(ncol),stat=alloc_stat) |
178 |
|
|
if (alloc_stat /= 0 ) then |
179 |
|
|
status = -1 |
180 |
|
|
return |
181 |
|
|
endif |
182 |
chuckv |
254 |
|
183 |
chuckv |
239 |
!! Gather idents into row and column idents |
184 |
chuckv |
254 |
|
185 |
chuckv |
239 |
call gather(ident,identRow,plan_row) |
186 |
|
|
call gather(ident,identCol,plan_col) |
187 |
chuckv |
254 |
|
188 |
chuckv |
239 |
|
189 |
|
|
!! Create row and col pointer lists |
190 |
chuckv |
254 |
|
191 |
|
|
call create_IdentPtrlst(identRow,ljListHead,identPtrListRow,thisStat) |
192 |
chuckv |
239 |
if (thisStat /= 0 ) then |
193 |
|
|
status = -1 |
194 |
|
|
return |
195 |
|
|
endif |
196 |
|
|
|
197 |
chuckv |
254 |
call create_IdentPtrlst(identCol,ljListHead,identPtrListColumn,thisStat) |
198 |
chuckv |
239 |
if (thisStat /= 0 ) then |
199 |
|
|
status = -1 |
200 |
|
|
return |
201 |
|
|
endif |
202 |
|
|
|
203 |
|
|
!! free temporary ident arrays |
204 |
|
|
deallocate(identCol) |
205 |
|
|
deallocate(identRow) |
206 |
|
|
|
207 |
|
|
|
208 |
|
|
!! Allocate neighbor lists for mpi simulations. |
209 |
|
|
if (.not. allocated(point)) then |
210 |
|
|
allocate(point(nrow),stat=alloc_stat) |
211 |
|
|
if (alloc_stat /=0) then |
212 |
|
|
status = -1 |
213 |
|
|
return |
214 |
|
|
endif |
215 |
|
|
|
216 |
|
|
allocate(list(ncol*listMultiplier),stat=alloc_stat) |
217 |
|
|
if (alloc_stat /=0) then |
218 |
|
|
status = -1 |
219 |
|
|
return |
220 |
|
|
endif |
221 |
|
|
else |
222 |
|
|
deallocate(list) |
223 |
|
|
deallocate(point) |
224 |
|
|
allocate(point(nrow),stat=alloc_stat) |
225 |
|
|
if (alloc_stat /=0) then |
226 |
|
|
status = -1 |
227 |
|
|
return |
228 |
|
|
endif |
229 |
|
|
|
230 |
|
|
allocate(list(ncol*listMultiplier),stat=alloc_stat) |
231 |
|
|
if (alloc_stat /=0) then |
232 |
|
|
status = -1 |
233 |
|
|
return |
234 |
|
|
endif |
235 |
|
|
endif |
236 |
|
|
|
237 |
chuckv |
230 |
#endif |
238 |
chuckv |
239 |
|
239 |
chuckv |
247 |
call createMixingList(thisStat) |
240 |
|
|
if (thisStat /= 0) then |
241 |
|
|
status = -1 |
242 |
|
|
return |
243 |
|
|
endif |
244 |
chuckv |
239 |
isljFFinit = .true. |
245 |
|
|
|
246 |
chuckv |
230 |
end subroutine init_ljFF |
247 |
|
|
|
248 |
chuckv |
239 |
|
249 |
|
|
|
250 |
|
|
|
251 |
|
|
|
252 |
chuckv |
230 |
|
253 |
chuckv |
247 |
subroutine createMixingList(status) |
254 |
|
|
integer :: listSize |
255 |
|
|
integer :: status |
256 |
chuckv |
230 |
integer :: i |
257 |
chuckv |
247 |
integer :: j |
258 |
chuckv |
230 |
|
259 |
chuckv |
247 |
integer :: outerCounter = 0 |
260 |
|
|
integer :: innerCounter = 0 |
261 |
|
|
type (lj_atype), pointer :: tmpPtrOuter => null() |
262 |
|
|
type (lj_atype), pointer :: tmpPtrInner => null() |
263 |
|
|
status = 0 |
264 |
chuckv |
230 |
|
265 |
chuckv |
247 |
listSize = getListLen(ljListHead) |
266 |
|
|
if (listSize == 0) then |
267 |
|
|
status = -1 |
268 |
|
|
return |
269 |
|
|
end if |
270 |
|
|
|
271 |
|
|
|
272 |
|
|
if (.not. associated(ljMixed)) then |
273 |
|
|
allocate(ljMixed(listSize,listSize)) |
274 |
chuckv |
230 |
else |
275 |
chuckv |
247 |
status = -1 |
276 |
|
|
return |
277 |
|
|
end if |
278 |
chuckv |
230 |
|
279 |
chuckv |
247 |
|
280 |
chuckv |
252 |
|
281 |
chuckv |
247 |
tmpPtrOuter => ljListHead |
282 |
|
|
tmpPtrInner => tmpPtrOuter%next |
283 |
|
|
do while (associated(tmpPtrOuter)) |
284 |
|
|
outerCounter = outerCounter + 1 |
285 |
|
|
! do self mixing rule |
286 |
|
|
ljMixed(outerCounter,outerCounter)%sigma = tmpPtrOuter%sigma |
287 |
|
|
|
288 |
|
|
ljMixed(outerCounter,outerCounter)%sigma2 = ljMixed(outerCounter,outerCounter)%sigma & |
289 |
|
|
* ljMixed(outerCounter,outerCounter)%sigma |
290 |
|
|
|
291 |
|
|
ljMixed(outerCounter,outerCounter)%sigma6 = ljMixed(outerCounter,outerCounter)%sigma2 & |
292 |
|
|
* ljMixed(outerCounter,outerCounter)%sigma2 & |
293 |
|
|
* ljMixed(outerCounter,outerCounter)%sigma2 |
294 |
|
|
|
295 |
|
|
ljMixed(outerCounter,outerCounter)%epsilon = tmpPtrOuter%epsilon |
296 |
|
|
|
297 |
|
|
innerCounter = outerCounter + 1 |
298 |
|
|
do while (associated(tmpPtrInner)) |
299 |
|
|
|
300 |
|
|
ljMixed(outerCounter,innerCounter)%sigma = & |
301 |
|
|
calcLJMix("sigma",tmpPtrOuter%sigma, & |
302 |
|
|
tmpPtrInner%sigma) |
303 |
|
|
|
304 |
|
|
ljMixed(outerCounter,innerCounter)%sigma2 = & |
305 |
|
|
ljMixed(outerCounter,innerCounter)%sigma & |
306 |
|
|
* ljMixed(outerCounter,innerCounter)%sigma |
307 |
|
|
|
308 |
|
|
ljMixed(outerCounter,innerCounter)%sigma6 = & |
309 |
|
|
ljMixed(outerCounter,innerCounter)%sigma2 & |
310 |
|
|
* ljMixed(outerCounter,innerCounter)%sigma2 & |
311 |
|
|
* ljMixed(outerCounter,innerCounter)%sigma2 |
312 |
|
|
|
313 |
|
|
ljMixed(outerCounter,innerCounter)%epsilon = & |
314 |
|
|
calcLJMix("epsilon",tmpPtrOuter%epsilon, & |
315 |
|
|
tmpPtrInner%epsilon) |
316 |
|
|
ljMixed(innerCounter,outerCounter)%sigma = ljMixed(outerCounter,innerCounter)%sigma |
317 |
|
|
ljMixed(innerCounter,outerCounter)%sigma2 = ljMixed(outerCounter,innerCounter)%sigma2 |
318 |
|
|
ljMixed(innerCounter,outerCounter)%sigma6 = ljMixed(outerCounter,innerCounter)%sigma6 |
319 |
|
|
ljMixed(innerCounter,outerCounter)%epsilon = ljMixed(outerCounter,innerCounter)%epsilon |
320 |
|
|
|
321 |
|
|
|
322 |
|
|
tmpPtrInner => tmpPtrInner%next |
323 |
|
|
innerCounter = innerCounter + 1 |
324 |
|
|
end do |
325 |
|
|
! advance pointers |
326 |
|
|
tmpPtrOuter => tmpPtrOuter%next |
327 |
|
|
if (associated(tmpPtrOuter)) then |
328 |
|
|
tmpPtrInner => tmpPtrOuter%next |
329 |
chuckv |
230 |
endif |
330 |
|
|
|
331 |
|
|
end do |
332 |
|
|
|
333 |
chuckv |
247 |
end subroutine createMixingList |
334 |
chuckv |
230 |
|
335 |
|
|
|
336 |
chuckv |
246 |
|
337 |
chuckv |
230 |
|
338 |
|
|
|
339 |
chuckv |
246 |
|
340 |
|
|
|
341 |
|
|
|
342 |
chuckv |
239 |
!! FORCE routine Calculates Lennard Jones forces. |
343 |
chuckv |
230 |
!-------------------------------------------------------------> |
344 |
|
|
subroutine do_lj_ff(q,f,potE,do_pot) |
345 |
chuckv |
246 |
!! Position array provided by C, dimensioned by getNlocal |
346 |
|
|
real ( kind = dp ), dimension(3,getNlocal()) :: q |
347 |
|
|
!! Force array provided by C, dimensioned by getNlocal |
348 |
|
|
real ( kind = dp ), dimension(3,getNlocal()) :: f |
349 |
chuckv |
222 |
real ( kind = dp ) :: potE |
350 |
mmeineke |
235 |
logical ( kind = 2) :: do_pot |
351 |
chuckv |
239 |
|
352 |
|
|
type(lj_atype), pointer :: ljAtype_i |
353 |
|
|
type(lj_atype), pointer :: ljAtype_j |
354 |
|
|
|
355 |
chuckv |
253 |
#ifdef IS_MPI |
356 |
chuckv |
254 |
real( kind = DP ), dimension(3,getNcol(plan_col)) :: efr |
357 |
chuckv |
230 |
real( kind = DP ) :: pot_local |
358 |
|
|
#else |
359 |
chuckv |
246 |
real( kind = DP ), dimension(3,getNlocal()) :: efr |
360 |
chuckv |
230 |
#endif |
361 |
|
|
|
362 |
chuckv |
247 |
!! Local arrays needed for MPI |
363 |
|
|
#ifdef IS_MPI |
364 |
chuckv |
254 |
real(kind = dp), dimension(3,getNrow(plan_row)) :: qRow |
365 |
|
|
real(kind = dp), dimension(3,getNcol(plan_col)) :: qCol |
366 |
chuckv |
247 |
|
367 |
chuckv |
254 |
real(kind = dp), dimension(3,getNrow(plan_row)) :: fRow |
368 |
|
|
real(kind = dp), dimension(3,getNcol(plan_col)) :: fCol |
369 |
|
|
real(kind = dp), dimension(3,getNlocal()) :: fMPITemp |
370 |
chuckv |
247 |
|
371 |
chuckv |
254 |
real(kind = dp), dimension(getNrow(plan_row)) :: eRow |
372 |
|
|
real(kind = dp), dimension(getNcol(plan_col)) :: eCol |
373 |
chuckv |
247 |
|
374 |
|
|
real(kind = dp), dimension(getNlocal()) :: eTemp |
375 |
|
|
#endif |
376 |
|
|
|
377 |
|
|
|
378 |
|
|
|
379 |
chuckv |
230 |
real( kind = DP ) :: pe |
380 |
chuckv |
246 |
logical :: update_nlist |
381 |
chuckv |
216 |
|
382 |
chuckv |
232 |
|
383 |
chuckv |
230 |
integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2 |
384 |
|
|
integer :: nlist |
385 |
|
|
integer :: j_start |
386 |
|
|
integer :: tag_i,tag_j |
387 |
|
|
real( kind = DP ) :: r, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp |
388 |
|
|
real( kind = DP ) :: rxi, ryi, rzi, rxij, ryij, rzij, rijsq |
389 |
chuckv |
246 |
real( kind = DP ) :: rlistsq, rcutsq,rlist,rcut |
390 |
chuckv |
222 |
|
391 |
chuckv |
246 |
! a rig that need to be fixed. |
392 |
|
|
#ifdef IS_MPI |
393 |
|
|
logical :: newtons_thrd = .true. |
394 |
|
|
real( kind = dp ) :: pe_local |
395 |
|
|
integer :: nlocal |
396 |
|
|
#endif |
397 |
chuckv |
230 |
integer :: nrow |
398 |
|
|
integer :: ncol |
399 |
chuckv |
246 |
integer :: natoms |
400 |
chuckv |
222 |
|
401 |
chuckv |
253 |
|
402 |
|
|
|
403 |
|
|
|
404 |
chuckv |
246 |
! Make sure we are properly initialized. |
405 |
chuckv |
239 |
if (.not. isljFFInit) then |
406 |
|
|
write(default_error,*) "ERROR: lj_FF has not been properly initialized" |
407 |
|
|
return |
408 |
|
|
endif |
409 |
chuckv |
246 |
#ifdef IS_MPI |
410 |
|
|
if (.not. isMPISimSet()) then |
411 |
|
|
write(default_error,*) "ERROR: mpiSimulation has not been properly initialized" |
412 |
|
|
return |
413 |
|
|
endif |
414 |
|
|
#endif |
415 |
chuckv |
222 |
|
416 |
chuckv |
246 |
!! initialize local variables |
417 |
|
|
natoms = getNlocal() |
418 |
|
|
call getRcut(rcut,rcut2=rcutsq) |
419 |
|
|
call getRlist(rlist,rlistsq) |
420 |
chuckv |
253 |
|
421 |
chuckv |
230 |
#ifndef IS_MPI |
422 |
|
|
nrow = natoms - 1 |
423 |
|
|
ncol = natoms |
424 |
|
|
#else |
425 |
chuckv |
239 |
nrow = getNrow(plan_row) |
426 |
|
|
ncol = getNcol(plan_col) |
427 |
chuckv |
246 |
nlocal = natoms |
428 |
chuckv |
230 |
j_start = 1 |
429 |
|
|
#endif |
430 |
chuckv |
222 |
|
431 |
chuckv |
239 |
|
432 |
chuckv |
246 |
!! See if we need to update neighbor lists |
433 |
|
|
call check(q,update_nlist) |
434 |
chuckv |
253 |
if (firstTime) then |
435 |
|
|
update_nlist = .true. |
436 |
|
|
firstTime = .false. |
437 |
|
|
endif |
438 |
chuckv |
239 |
|
439 |
|
|
!--------------WARNING........................... |
440 |
|
|
! Zero variables, NOTE:::: Forces are zeroed in C |
441 |
|
|
! Zeroing them here could delete previously computed |
442 |
|
|
! Forces. |
443 |
|
|
!------------------------------------------------ |
444 |
chuckv |
230 |
#ifndef IS_MPI |
445 |
chuckv |
246 |
! nloops = nloops + 1 |
446 |
|
|
pe = 0.0E0_DP |
447 |
|
|
|
448 |
chuckv |
230 |
#else |
449 |
chuckv |
253 |
fRow = 0.0E0_DP |
450 |
|
|
fCol = 0.0E0_DP |
451 |
chuckv |
230 |
|
452 |
chuckv |
246 |
pe_local = 0.0E0_DP |
453 |
chuckv |
230 |
|
454 |
chuckv |
246 |
eRow = 0.0E0_DP |
455 |
|
|
eCol = 0.0E0_DP |
456 |
|
|
eTemp = 0.0E0_DP |
457 |
chuckv |
230 |
#endif |
458 |
|
|
efr = 0.0E0_DP |
459 |
|
|
|
460 |
|
|
! communicate MPI positions |
461 |
chuckv |
253 |
#ifdef IS_MPI |
462 |
|
|
call gather(q,qRow,plan_row3d) |
463 |
|
|
call gather(q,qCol,plan_col3d) |
464 |
chuckv |
230 |
#endif |
465 |
|
|
|
466 |
|
|
|
467 |
|
|
if (update_nlist) then |
468 |
|
|
|
469 |
|
|
! save current configuration, contruct neighbor list, |
470 |
|
|
! and calculate forces |
471 |
chuckv |
246 |
call save_nlist(q) |
472 |
chuckv |
230 |
|
473 |
|
|
nlist = 0 |
474 |
|
|
|
475 |
chuckv |
253 |
|
476 |
chuckv |
230 |
|
477 |
|
|
do i = 1, nrow |
478 |
|
|
point(i) = nlist + 1 |
479 |
chuckv |
253 |
#ifdef IS_MPI |
480 |
chuckv |
239 |
ljAtype_i => identPtrListRow(i)%this |
481 |
|
|
tag_i = tagRow(i) |
482 |
|
|
rxi = qRow(1,i) |
483 |
|
|
ryi = qRow(2,i) |
484 |
|
|
rzi = qRow(3,i) |
485 |
chuckv |
230 |
#else |
486 |
chuckv |
239 |
ljAtype_i => identPtrList(i)%this |
487 |
chuckv |
230 |
j_start = i + 1 |
488 |
|
|
rxi = q(1,i) |
489 |
|
|
ryi = q(2,i) |
490 |
|
|
rzi = q(3,i) |
491 |
|
|
#endif |
492 |
|
|
|
493 |
|
|
inner: do j = j_start, ncol |
494 |
chuckv |
253 |
#ifdef IS_MPI |
495 |
chuckv |
239 |
! Assign identity pointers and tags |
496 |
|
|
ljAtype_j => identPtrListColumn(j)%this |
497 |
chuckv |
254 |
tag_j = tagColumn(j) |
498 |
chuckv |
230 |
if (newtons_thrd) then |
499 |
|
|
if (tag_i <= tag_j) then |
500 |
|
|
if (mod(tag_i + tag_j,2) == 0) cycle inner |
501 |
|
|
else |
502 |
|
|
if (mod(tag_i + tag_j,2) == 1) cycle inner |
503 |
|
|
endif |
504 |
|
|
endif |
505 |
|
|
|
506 |
chuckv |
239 |
rxij = wrap(rxi - qCol(1,j), 1) |
507 |
|
|
ryij = wrap(ryi - qCol(2,j), 2) |
508 |
|
|
rzij = wrap(rzi - qCol(3,j), 3) |
509 |
chuckv |
230 |
#else |
510 |
chuckv |
239 |
ljAtype_j => identPtrList(j)%this |
511 |
chuckv |
230 |
rxij = wrap(rxi - q(1,j), 1) |
512 |
|
|
ryij = wrap(ryi - q(2,j), 2) |
513 |
|
|
rzij = wrap(rzi - q(3,j), 3) |
514 |
|
|
|
515 |
|
|
#endif |
516 |
|
|
rijsq = rxij*rxij + ryij*ryij + rzij*rzij |
517 |
|
|
|
518 |
chuckv |
253 |
#ifdef IS_MPI |
519 |
chuckv |
246 |
if (rijsq <= rlistsq .AND. & |
520 |
chuckv |
230 |
tag_j /= tag_i) then |
521 |
|
|
#else |
522 |
chuckv |
253 |
|
523 |
chuckv |
246 |
if (rijsq < rlistsq) then |
524 |
chuckv |
230 |
#endif |
525 |
|
|
|
526 |
|
|
nlist = nlist + 1 |
527 |
|
|
if (nlist > size(list)) then |
528 |
chuckv |
246 |
!! "Change how nlist size is done" |
529 |
chuckv |
239 |
write(DEFAULT_ERROR,*) "ERROR: nlist > list size" |
530 |
chuckv |
230 |
endif |
531 |
|
|
list(nlist) = j |
532 |
|
|
|
533 |
chuckv |
253 |
|
534 |
chuckv |
230 |
if (rijsq < rcutsq) then |
535 |
|
|
|
536 |
|
|
r = dsqrt(rijsq) |
537 |
|
|
|
538 |
chuckv |
239 |
call getLJPot(r,pot,dudr,ljAtype_i,ljAtype_j) |
539 |
chuckv |
230 |
|
540 |
chuckv |
253 |
#ifdef IS_MPI |
541 |
chuckv |
246 |
eRow(i) = eRow(i) + pot*0.5 |
542 |
|
|
eCol(i) = eCol(i) + pot*0.5 |
543 |
chuckv |
230 |
#else |
544 |
chuckv |
253 |
pe = pe + pot |
545 |
chuckv |
230 |
#endif |
546 |
|
|
|
547 |
|
|
efr(1,j) = -rxij |
548 |
|
|
efr(2,j) = -ryij |
549 |
|
|
efr(3,j) = -rzij |
550 |
|
|
|
551 |
|
|
do dim = 1, 3 |
552 |
|
|
|
553 |
|
|
|
554 |
|
|
drdx1 = efr(dim,j) / r |
555 |
|
|
ftmp = dudr * drdx1 |
556 |
|
|
|
557 |
|
|
|
558 |
chuckv |
253 |
#ifdef IS_MPI |
559 |
chuckv |
239 |
fCol(dim,j) = fCol(dim,j) - ftmp |
560 |
|
|
fRow(dim,i) = fRow(dim,i) + ftmp |
561 |
chuckv |
230 |
#else |
562 |
|
|
|
563 |
|
|
f(dim,j) = f(dim,j) - ftmp |
564 |
|
|
f(dim,i) = f(dim,i) + ftmp |
565 |
|
|
|
566 |
|
|
#endif |
567 |
|
|
enddo |
568 |
|
|
endif |
569 |
|
|
endif |
570 |
|
|
enddo inner |
571 |
|
|
enddo |
572 |
|
|
|
573 |
chuckv |
253 |
#ifdef IS_MPI |
574 |
chuckv |
230 |
point(nrow + 1) = nlist + 1 |
575 |
|
|
#else |
576 |
|
|
point(natoms) = nlist + 1 |
577 |
|
|
#endif |
578 |
|
|
|
579 |
|
|
else |
580 |
|
|
|
581 |
|
|
! use the list to find the neighbors |
582 |
|
|
do i = 1, nrow |
583 |
|
|
JBEG = POINT(i) |
584 |
|
|
JEND = POINT(i+1) - 1 |
585 |
|
|
! check thiat molecule i has neighbors |
586 |
|
|
if (jbeg .le. jend) then |
587 |
chuckv |
253 |
#ifdef IS_MPI |
588 |
chuckv |
239 |
ljAtype_i => identPtrListRow(i)%this |
589 |
|
|
rxi = qRow(1,i) |
590 |
|
|
ryi = qRow(2,i) |
591 |
|
|
rzi = qRow(3,i) |
592 |
chuckv |
230 |
#else |
593 |
chuckv |
239 |
ljAtype_i => identPtrList(i)%this |
594 |
chuckv |
230 |
rxi = q(1,i) |
595 |
|
|
ryi = q(2,i) |
596 |
|
|
rzi = q(3,i) |
597 |
|
|
#endif |
598 |
|
|
do jnab = jbeg, jend |
599 |
|
|
j = list(jnab) |
600 |
chuckv |
253 |
#ifdef IS_MPI |
601 |
chuckv |
239 |
ljAtype_j = identPtrListColumn(j)%this |
602 |
chuckv |
254 |
rxij = wrap(rxi - qCol(1,j), 1) |
603 |
|
|
ryij = wrap(ryi - qCol(2,j), 2) |
604 |
|
|
rzij = wrap(rzi - qCol(3,j), 3) |
605 |
chuckv |
230 |
#else |
606 |
chuckv |
239 |
ljAtype_j = identPtrList(j)%this |
607 |
|
|
rxij = wrap(rxi - q(1,j), 1) |
608 |
|
|
ryij = wrap(ryi - q(2,j), 2) |
609 |
|
|
rzij = wrap(rzi - q(3,j), 3) |
610 |
chuckv |
230 |
#endif |
611 |
|
|
rijsq = rxij*rxij + ryij*ryij + rzij*rzij |
612 |
|
|
|
613 |
|
|
if (rijsq < rcutsq) then |
614 |
|
|
|
615 |
|
|
r = dsqrt(rijsq) |
616 |
|
|
|
617 |
chuckv |
239 |
call getLJPot(r,pot,dudr,ljAtype_i,ljAtype_j) |
618 |
chuckv |
253 |
#ifdef IS_MPI |
619 |
chuckv |
246 |
eRow(i) = eRow(i) + pot*0.5 |
620 |
|
|
eCol(i) = eCol(i) + pot*0.5 |
621 |
chuckv |
230 |
#else |
622 |
chuckv |
246 |
pe = pe + pot |
623 |
chuckv |
230 |
#endif |
624 |
|
|
|
625 |
|
|
|
626 |
|
|
efr(1,j) = -rxij |
627 |
|
|
efr(2,j) = -ryij |
628 |
|
|
efr(3,j) = -rzij |
629 |
|
|
|
630 |
|
|
do dim = 1, 3 |
631 |
|
|
|
632 |
|
|
drdx1 = efr(dim,j) / r |
633 |
|
|
ftmp = dudr * drdx1 |
634 |
chuckv |
253 |
#ifdef IS_MPI |
635 |
chuckv |
239 |
fCol(dim,j) = fCol(dim,j) - ftmp |
636 |
|
|
fRow(dim,i) = fRow(dim,i) + ftmp |
637 |
chuckv |
230 |
#else |
638 |
|
|
f(dim,j) = f(dim,j) - ftmp |
639 |
|
|
f(dim,i) = f(dim,i) + ftmp |
640 |
|
|
#endif |
641 |
|
|
enddo |
642 |
|
|
endif |
643 |
|
|
enddo |
644 |
|
|
endif |
645 |
|
|
enddo |
646 |
|
|
endif |
647 |
|
|
|
648 |
|
|
|
649 |
|
|
|
650 |
chuckv |
253 |
#ifdef IS_MPI |
651 |
chuckv |
230 |
!!distribute forces |
652 |
chuckv |
253 |
call scatter(fRow,f,plan_row3d) |
653 |
chuckv |
230 |
|
654 |
chuckv |
253 |
call scatter(fCol,fMPITemp,plan_col3d) |
655 |
chuckv |
246 |
|
656 |
chuckv |
230 |
do i = 1,nlocal |
657 |
chuckv |
253 |
f(1:3,i) = f(1:3,i) + fMPITemp(1:3,i) |
658 |
chuckv |
230 |
end do |
659 |
|
|
|
660 |
|
|
|
661 |
|
|
|
662 |
|
|
if (do_pot) then |
663 |
|
|
! scatter/gather pot_row into the members of my column |
664 |
chuckv |
253 |
call scatter(eRow,eTemp,plan_row) |
665 |
chuckv |
230 |
|
666 |
|
|
! scatter/gather pot_local into all other procs |
667 |
|
|
! add resultant to get total pot |
668 |
|
|
do i = 1, nlocal |
669 |
chuckv |
253 |
pe_local = pe_local + eTemp(i) |
670 |
chuckv |
230 |
enddo |
671 |
|
|
if (newtons_thrd) then |
672 |
chuckv |
253 |
eTemp = 0.0E0_DP |
673 |
|
|
call scatter(eCol,eTemp,plan_col) |
674 |
chuckv |
230 |
do i = 1, nlocal |
675 |
chuckv |
253 |
pe_local = pe_local + eTemp(i) |
676 |
chuckv |
230 |
enddo |
677 |
|
|
endif |
678 |
|
|
endif |
679 |
chuckv |
246 |
potE = pe_local |
680 |
chuckv |
230 |
#endif |
681 |
|
|
|
682 |
chuckv |
246 |
potE = pe |
683 |
chuckv |
230 |
|
684 |
chuckv |
253 |
|
685 |
chuckv |
230 |
|
686 |
chuckv |
253 |
|
687 |
chuckv |
222 |
end subroutine do_lj_ff |
688 |
|
|
|
689 |
chuckv |
239 |
!! Calculates the potential between two lj particles based on two lj_atype pointers, optionally returns second |
690 |
chuckv |
230 |
!! derivatives. |
691 |
|
|
subroutine getLjPot(r,pot,dudr,atype1,atype2,d2,status) |
692 |
|
|
! arguments |
693 |
|
|
!! Length of vector between particles |
694 |
|
|
real( kind = dp ), intent(in) :: r |
695 |
|
|
!! Potential Energy |
696 |
|
|
real( kind = dp ), intent(out) :: pot |
697 |
|
|
!! Derivatve wrt postion |
698 |
|
|
real( kind = dp ), intent(out) :: dudr |
699 |
|
|
!! Second Derivative, optional, used mainly for normal mode calculations. |
700 |
|
|
real( kind = dp ), intent(out), optional :: d2 |
701 |
|
|
|
702 |
chuckv |
246 |
type (lj_atype), pointer :: atype1 |
703 |
|
|
type (lj_atype), pointer :: atype2 |
704 |
chuckv |
222 |
|
705 |
chuckv |
230 |
integer, intent(out), optional :: status |
706 |
chuckv |
222 |
|
707 |
chuckv |
230 |
! local Variables |
708 |
|
|
real( kind = dp ) :: sigma |
709 |
|
|
real( kind = dp ) :: sigma2 |
710 |
|
|
real( kind = dp ) :: sigma6 |
711 |
chuckv |
246 |
real( kind = dp ) :: epsilon |
712 |
chuckv |
230 |
|
713 |
|
|
real( kind = dp ) :: rcut |
714 |
|
|
real( kind = dp ) :: rcut2 |
715 |
|
|
real( kind = dp ) :: rcut6 |
716 |
|
|
real( kind = dp ) :: r2 |
717 |
|
|
real( kind = dp ) :: r6 |
718 |
|
|
|
719 |
|
|
real( kind = dp ) :: t6 |
720 |
|
|
real( kind = dp ) :: t12 |
721 |
|
|
real( kind = dp ) :: tp6 |
722 |
|
|
real( kind = dp ) :: tp12 |
723 |
|
|
real( kind = dp ) :: delta |
724 |
|
|
|
725 |
|
|
logical :: doSec = .false. |
726 |
|
|
|
727 |
|
|
integer :: errorStat |
728 |
|
|
|
729 |
|
|
!! Optional Argument Checking |
730 |
|
|
! Check to see if we need to do second derivatives |
731 |
|
|
|
732 |
|
|
if (present(d2)) doSec = .true. |
733 |
|
|
if (present(status)) status = 0 |
734 |
|
|
|
735 |
|
|
! Look up the correct parameters in the mixing matrix |
736 |
chuckv |
246 |
sigma = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma |
737 |
|
|
sigma2 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma2 |
738 |
|
|
sigma6 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma6 |
739 |
|
|
epsilon = ljMixed(atype1%atype_ident,atype2%atype_ident)%epsilon |
740 |
chuckv |
230 |
|
741 |
|
|
|
742 |
chuckv |
253 |
|
743 |
|
|
|
744 |
chuckv |
230 |
call getRcut(rcut,rcut2=rcut2,rcut6=rcut6,status=errorStat) |
745 |
|
|
|
746 |
|
|
r2 = r * r |
747 |
|
|
r6 = r2 * r2 * r2 |
748 |
|
|
|
749 |
|
|
t6 = sigma6/ r6 |
750 |
|
|
t12 = t6 * t6 |
751 |
|
|
|
752 |
|
|
|
753 |
|
|
|
754 |
|
|
tp6 = sigma6 / rcut6 |
755 |
|
|
tp12 = tp6*tp6 |
756 |
|
|
|
757 |
|
|
delta = -4.0E0_DP*epsilon * (tp12 - tp6) |
758 |
|
|
|
759 |
|
|
if (r.le.rcut) then |
760 |
chuckv |
246 |
pot = 4.0E0_DP * epsilon * (t12 - t6) + delta |
761 |
chuckv |
230 |
dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / r |
762 |
|
|
if(doSec) d2 = 24.0E0_DP * epsilon * (26.0E0_DP*t12 - 7.0E0_DP*t6)/r/r |
763 |
|
|
else |
764 |
chuckv |
246 |
pot = 0.0E0_DP |
765 |
chuckv |
230 |
dudr = 0.0E0_DP |
766 |
|
|
if(doSec) d2 = 0.0E0_DP |
767 |
|
|
endif |
768 |
|
|
|
769 |
|
|
return |
770 |
|
|
|
771 |
|
|
|
772 |
|
|
|
773 |
|
|
end subroutine getLjPot |
774 |
|
|
|
775 |
|
|
|
776 |
chuckv |
239 |
!! Calculates the mixing for sigma or epslon based on character string for initialzition of mixing array. |
777 |
chuckv |
230 |
function calcLJMix(thisParam,param1,param2,status) result(myMixParam) |
778 |
|
|
character(len=*) :: thisParam |
779 |
|
|
real(kind = dp) :: param1 |
780 |
|
|
real(kind = dp) :: param2 |
781 |
|
|
real(kind = dp ) :: myMixParam |
782 |
|
|
integer, optional :: status |
783 |
|
|
|
784 |
|
|
|
785 |
|
|
myMixParam = 0.0_dp |
786 |
|
|
|
787 |
|
|
if (present(status)) status = 0 |
788 |
|
|
|
789 |
|
|
select case (thisParam) |
790 |
|
|
|
791 |
|
|
case ("sigma") |
792 |
|
|
myMixParam = 0.5_dp * (param1 + param2) |
793 |
chuckv |
246 |
case ("epsilon") |
794 |
chuckv |
230 |
myMixParam = sqrt(param1 * param2) |
795 |
|
|
case default |
796 |
|
|
status = -1 |
797 |
|
|
end select |
798 |
|
|
|
799 |
|
|
end function calcLJMix |
800 |
|
|
|
801 |
|
|
|
802 |
|
|
|
803 |
chuckv |
215 |
end module lj_ff |