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