1 |
mmeineke |
377 |
!! do_Forces.F90 |
2 |
|
|
!! module do_Forces |
3 |
|
|
!! Calculates Long Range forces. |
4 |
|
|
|
5 |
|
|
!! @author Charles F. Vardeman II |
6 |
|
|
!! @author Matthew Meineke |
7 |
mmeineke |
634 |
!! @version $Id: do_Forces.F90,v 1.23 2003-07-17 19:38:23 mmeineke Exp $, $Date: 2003-07-17 19:38:23 $, $Name: not supported by cvs2svn $, $Revision: 1.23 $ |
8 |
mmeineke |
377 |
|
9 |
|
|
module do_Forces |
10 |
|
|
use force_globals |
11 |
|
|
use simulation |
12 |
|
|
use definitions |
13 |
|
|
use atype_module |
14 |
|
|
use neighborLists |
15 |
|
|
use lj |
16 |
|
|
use sticky_pair |
17 |
|
|
use dipole_dipole |
18 |
|
|
use reaction_field |
19 |
|
|
use gb_pair |
20 |
mmeineke |
626 |
use vector_class |
21 |
mmeineke |
377 |
#ifdef IS_MPI |
22 |
|
|
use mpiSimulation |
23 |
|
|
#endif |
24 |
|
|
|
25 |
|
|
implicit none |
26 |
|
|
PRIVATE |
27 |
|
|
|
28 |
|
|
#define __FORTRAN90 |
29 |
|
|
#include "fForceField.h" |
30 |
|
|
|
31 |
mmeineke |
626 |
logical, save :: do_forces_initialized = .false., haveRlist = .false. |
32 |
|
|
logical, save :: havePolicies = .false. |
33 |
mmeineke |
377 |
logical, save :: FF_uses_LJ |
34 |
|
|
logical, save :: FF_uses_sticky |
35 |
|
|
logical, save :: FF_uses_dipoles |
36 |
|
|
logical, save :: FF_uses_RF |
37 |
|
|
logical, save :: FF_uses_GB |
38 |
|
|
logical, save :: FF_uses_EAM |
39 |
|
|
|
40 |
mmeineke |
626 |
real(kind=dp), save :: rlist, rlistsq |
41 |
|
|
|
42 |
mmeineke |
377 |
public :: init_FF |
43 |
|
|
public :: do_force_loop |
44 |
mmeineke |
626 |
public :: setRlistDF |
45 |
mmeineke |
377 |
|
46 |
|
|
contains |
47 |
|
|
|
48 |
mmeineke |
626 |
subroutine setRlistDF( this_rlist ) |
49 |
|
|
|
50 |
|
|
real(kind=dp) :: this_rlist |
51 |
|
|
|
52 |
|
|
rlist = this_rlist |
53 |
|
|
rlistsq = rlist * rlist |
54 |
|
|
|
55 |
|
|
haveRlist = .true. |
56 |
|
|
if( havePolicies ) do_forces_initialized = .true. |
57 |
|
|
|
58 |
|
|
end subroutine setRlistDF |
59 |
|
|
|
60 |
mmeineke |
377 |
subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat) |
61 |
|
|
|
62 |
|
|
integer, intent(in) :: LJMIXPOLICY |
63 |
|
|
logical, intent(in) :: use_RF_c |
64 |
|
|
|
65 |
|
|
integer, intent(out) :: thisStat |
66 |
|
|
integer :: my_status, nMatches |
67 |
|
|
integer, pointer :: MatchList(:) => null() |
68 |
|
|
real(kind=dp) :: rcut, rrf, rt, dielect |
69 |
|
|
|
70 |
|
|
!! assume things are copacetic, unless they aren't |
71 |
|
|
thisStat = 0 |
72 |
|
|
|
73 |
|
|
!! Fortran's version of a cast: |
74 |
|
|
FF_uses_RF = use_RF_c |
75 |
|
|
|
76 |
|
|
!! init_FF is called *after* all of the atom types have been |
77 |
|
|
!! defined in atype_module using the new_atype subroutine. |
78 |
|
|
!! |
79 |
|
|
!! this will scan through the known atypes and figure out what |
80 |
|
|
!! interactions are used by the force field. |
81 |
|
|
|
82 |
|
|
FF_uses_LJ = .false. |
83 |
|
|
FF_uses_sticky = .false. |
84 |
|
|
FF_uses_dipoles = .false. |
85 |
|
|
FF_uses_GB = .false. |
86 |
|
|
FF_uses_EAM = .false. |
87 |
|
|
|
88 |
|
|
call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList) |
89 |
|
|
if (nMatches .gt. 0) FF_uses_LJ = .true. |
90 |
|
|
|
91 |
|
|
call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList) |
92 |
|
|
if (nMatches .gt. 0) FF_uses_dipoles = .true. |
93 |
|
|
|
94 |
|
|
call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, & |
95 |
|
|
MatchList) |
96 |
|
|
if (nMatches .gt. 0) FF_uses_Sticky = .true. |
97 |
|
|
|
98 |
|
|
call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList) |
99 |
|
|
if (nMatches .gt. 0) FF_uses_GB = .true. |
100 |
|
|
|
101 |
|
|
call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList) |
102 |
|
|
if (nMatches .gt. 0) FF_uses_EAM = .true. |
103 |
|
|
|
104 |
|
|
!! check to make sure the FF_uses_RF setting makes sense |
105 |
|
|
|
106 |
|
|
if (FF_uses_dipoles) then |
107 |
|
|
if (FF_uses_RF) then |
108 |
|
|
dielect = getDielect() |
109 |
mmeineke |
626 |
call initialize_rf(dielect) |
110 |
mmeineke |
377 |
endif |
111 |
|
|
else |
112 |
|
|
if (FF_uses_RF) then |
113 |
|
|
write(default_error,*) 'Using Reaction Field with no dipoles? Huh?' |
114 |
|
|
thisStat = -1 |
115 |
|
|
return |
116 |
|
|
endif |
117 |
mmeineke |
626 |
endif |
118 |
mmeineke |
377 |
|
119 |
|
|
if (FF_uses_LJ) then |
120 |
|
|
|
121 |
|
|
select case (LJMIXPOLICY) |
122 |
|
|
case (LB_MIXING_RULE) |
123 |
mmeineke |
626 |
call init_lj_FF(LB_MIXING_RULE, my_status) |
124 |
mmeineke |
377 |
case (EXPLICIT_MIXING_RULE) |
125 |
mmeineke |
626 |
call init_lj_FF(EXPLICIT_MIXING_RULE, my_status) |
126 |
mmeineke |
377 |
case default |
127 |
|
|
write(default_error,*) 'unknown LJ Mixing Policy!' |
128 |
|
|
thisStat = -1 |
129 |
|
|
return |
130 |
|
|
end select |
131 |
|
|
if (my_status /= 0) then |
132 |
|
|
thisStat = -1 |
133 |
|
|
return |
134 |
|
|
end if |
135 |
|
|
endif |
136 |
|
|
|
137 |
|
|
if (FF_uses_sticky) then |
138 |
|
|
call check_sticky_FF(my_status) |
139 |
|
|
if (my_status /= 0) then |
140 |
|
|
thisStat = -1 |
141 |
|
|
return |
142 |
|
|
end if |
143 |
|
|
endif |
144 |
|
|
|
145 |
|
|
if (FF_uses_GB) then |
146 |
|
|
call check_gb_pair_FF(my_status) |
147 |
|
|
if (my_status .ne. 0) then |
148 |
|
|
thisStat = -1 |
149 |
|
|
return |
150 |
|
|
endif |
151 |
|
|
endif |
152 |
|
|
|
153 |
|
|
if (FF_uses_GB .and. FF_uses_LJ) then |
154 |
|
|
endif |
155 |
chuckv |
480 |
if (.not. do_forces_initialized) then |
156 |
|
|
!! Create neighbor lists |
157 |
|
|
call expandNeighborList(getNlocal(), my_status) |
158 |
|
|
if (my_Status /= 0) then |
159 |
|
|
write(default_error,*) "SimSetup: ExpandNeighborList returned error." |
160 |
|
|
thisStat = -1 |
161 |
|
|
return |
162 |
|
|
endif |
163 |
|
|
endif |
164 |
mmeineke |
377 |
|
165 |
mmeineke |
626 |
havePolicies = .true. |
166 |
|
|
if( haveRlist ) do_forces_initialized = .true. |
167 |
|
|
|
168 |
mmeineke |
377 |
end subroutine init_FF |
169 |
|
|
|
170 |
|
|
|
171 |
|
|
!! Does force loop over i,j pairs. Calls do_pair to calculates forces. |
172 |
|
|
!-------------------------------------------------------------> |
173 |
|
|
subroutine do_force_loop(q, A, u_l, f, t, tau, pot, do_pot_c, do_stress_c, & |
174 |
|
|
error) |
175 |
|
|
!! Position array provided by C, dimensioned by getNlocal |
176 |
|
|
real ( kind = dp ), dimension(3,getNlocal()) :: q |
177 |
|
|
!! Rotation Matrix for each long range particle in simulation. |
178 |
|
|
real( kind = dp), dimension(9,getNlocal()) :: A |
179 |
|
|
!! Unit vectors for dipoles (lab frame) |
180 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
181 |
|
|
!! Force array provided by C, dimensioned by getNlocal |
182 |
|
|
real ( kind = dp ), dimension(3,getNlocal()) :: f |
183 |
|
|
!! Torsion array provided by C, dimensioned by getNlocal |
184 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: t |
185 |
|
|
!! Stress Tensor |
186 |
|
|
real( kind = dp), dimension(9) :: tau |
187 |
|
|
real ( kind = dp ) :: pot |
188 |
|
|
logical ( kind = 2) :: do_pot_c, do_stress_c |
189 |
|
|
logical :: do_pot |
190 |
|
|
logical :: do_stress |
191 |
chuckv |
439 |
#ifdef IS_MPI |
192 |
chuckv |
441 |
real( kind = DP ) :: pot_local |
193 |
mmeineke |
377 |
integer :: nrow |
194 |
|
|
integer :: ncol |
195 |
|
|
#endif |
196 |
|
|
integer :: nlocal |
197 |
|
|
integer :: natoms |
198 |
|
|
logical :: update_nlist |
199 |
|
|
integer :: i, j, jbeg, jend, jnab |
200 |
|
|
integer :: nlist |
201 |
mmeineke |
626 |
real( kind = DP ) :: rijsq |
202 |
mmeineke |
377 |
real(kind=dp),dimension(3) :: d |
203 |
|
|
real(kind=dp) :: rfpot, mu_i, virial |
204 |
|
|
integer :: me_i |
205 |
|
|
logical :: is_dp_i |
206 |
|
|
integer :: neighborListSize |
207 |
|
|
integer :: listerror, error |
208 |
|
|
integer :: localError |
209 |
mmeineke |
626 |
|
210 |
|
|
real(kind=dp) :: listSkin = 1.0 |
211 |
mmeineke |
619 |
|
212 |
mmeineke |
377 |
|
213 |
|
|
!! initialize local variables |
214 |
|
|
|
215 |
|
|
#ifdef IS_MPI |
216 |
chuckv |
441 |
pot_local = 0.0_dp |
217 |
mmeineke |
377 |
nlocal = getNlocal() |
218 |
|
|
nrow = getNrow(plan_row) |
219 |
|
|
ncol = getNcol(plan_col) |
220 |
|
|
#else |
221 |
|
|
nlocal = getNlocal() |
222 |
|
|
natoms = nlocal |
223 |
|
|
#endif |
224 |
chuckv |
441 |
|
225 |
mmeineke |
377 |
call check_initialization(localError) |
226 |
|
|
if ( localError .ne. 0 ) then |
227 |
|
|
error = -1 |
228 |
|
|
return |
229 |
|
|
end if |
230 |
|
|
call zero_work_arrays() |
231 |
|
|
|
232 |
|
|
do_pot = do_pot_c |
233 |
|
|
do_stress = do_stress_c |
234 |
|
|
|
235 |
|
|
! Gather all information needed by all force loops: |
236 |
|
|
|
237 |
|
|
#ifdef IS_MPI |
238 |
|
|
|
239 |
|
|
call gather(q,q_Row,plan_row3d) |
240 |
|
|
call gather(q,q_Col,plan_col3d) |
241 |
|
|
|
242 |
|
|
if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then |
243 |
|
|
call gather(u_l,u_l_Row,plan_row3d) |
244 |
|
|
call gather(u_l,u_l_Col,plan_col3d) |
245 |
|
|
|
246 |
|
|
call gather(A,A_Row,plan_row_rotation) |
247 |
|
|
call gather(A,A_Col,plan_col_rotation) |
248 |
|
|
endif |
249 |
|
|
|
250 |
|
|
#endif |
251 |
|
|
|
252 |
|
|
if (FF_RequiresPrepairCalc() .and. SimRequiresPrepairCalc()) then |
253 |
|
|
!! See if we need to update neighbor lists |
254 |
mmeineke |
626 |
call checkNeighborList(nlocal, q, listSkin, update_nlist) |
255 |
mmeineke |
377 |
!! if_mpi_gather_stuff_for_prepair |
256 |
|
|
!! do_prepair_loop_if_needed |
257 |
|
|
!! if_mpi_scatter_stuff_from_prepair |
258 |
|
|
!! if_mpi_gather_stuff_from_prepair_to_main_loop |
259 |
|
|
else |
260 |
|
|
!! See if we need to update neighbor lists |
261 |
mmeineke |
626 |
call checkNeighborList(nlocal, q, listSkin, update_nlist) |
262 |
mmeineke |
377 |
endif |
263 |
|
|
|
264 |
|
|
#ifdef IS_MPI |
265 |
|
|
|
266 |
|
|
if (update_nlist) then |
267 |
|
|
|
268 |
|
|
!! save current configuration, construct neighbor list, |
269 |
|
|
!! and calculate forces |
270 |
mmeineke |
459 |
call saveNeighborList(nlocal, q) |
271 |
mmeineke |
377 |
|
272 |
|
|
neighborListSize = size(list) |
273 |
|
|
nlist = 0 |
274 |
|
|
|
275 |
|
|
do i = 1, nrow |
276 |
|
|
point(i) = nlist + 1 |
277 |
|
|
|
278 |
|
|
inner: do j = 1, ncol |
279 |
|
|
|
280 |
|
|
if (skipThisPair(i,j)) cycle inner |
281 |
|
|
|
282 |
|
|
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
283 |
|
|
|
284 |
mmeineke |
626 |
if (rijsq < rlistsq) then |
285 |
mmeineke |
377 |
|
286 |
|
|
nlist = nlist + 1 |
287 |
|
|
|
288 |
|
|
if (nlist > neighborListSize) then |
289 |
|
|
call expandNeighborList(nlocal, listerror) |
290 |
|
|
if (listerror /= 0) then |
291 |
|
|
error = -1 |
292 |
|
|
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
293 |
|
|
return |
294 |
|
|
end if |
295 |
|
|
neighborListSize = size(list) |
296 |
|
|
endif |
297 |
|
|
|
298 |
|
|
list(nlist) = j |
299 |
|
|
|
300 |
mmeineke |
626 |
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
301 |
|
|
u_l, A, f, t, pot_local) |
302 |
|
|
|
303 |
mmeineke |
377 |
endif |
304 |
|
|
enddo inner |
305 |
|
|
enddo |
306 |
|
|
|
307 |
|
|
point(nrow + 1) = nlist + 1 |
308 |
|
|
|
309 |
|
|
else !! (of update_check) |
310 |
|
|
|
311 |
|
|
! use the list to find the neighbors |
312 |
|
|
do i = 1, nrow |
313 |
|
|
JBEG = POINT(i) |
314 |
|
|
JEND = POINT(i+1) - 1 |
315 |
|
|
! check thiat molecule i has neighbors |
316 |
|
|
if (jbeg .le. jend) then |
317 |
|
|
|
318 |
|
|
do jnab = jbeg, jend |
319 |
|
|
j = list(jnab) |
320 |
|
|
|
321 |
|
|
call get_interatomic_vector(q_Row(:,i), q_Col(:,j), d, rijsq) |
322 |
|
|
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
323 |
chuckv |
441 |
u_l, A, f, t, pot_local) |
324 |
mmeineke |
377 |
|
325 |
|
|
enddo |
326 |
|
|
endif |
327 |
|
|
enddo |
328 |
|
|
endif |
329 |
|
|
|
330 |
|
|
#else |
331 |
|
|
|
332 |
|
|
if (update_nlist) then |
333 |
|
|
|
334 |
|
|
! save current configuration, contruct neighbor list, |
335 |
|
|
! and calculate forces |
336 |
mmeineke |
459 |
call saveNeighborList(natoms, q) |
337 |
mmeineke |
377 |
|
338 |
|
|
neighborListSize = size(list) |
339 |
|
|
|
340 |
|
|
nlist = 0 |
341 |
|
|
|
342 |
|
|
do i = 1, natoms-1 |
343 |
|
|
point(i) = nlist + 1 |
344 |
|
|
|
345 |
|
|
inner: do j = i+1, natoms |
346 |
|
|
|
347 |
chuckv |
388 |
if (skipThisPair(i,j)) cycle inner |
348 |
|
|
|
349 |
mmeineke |
377 |
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
350 |
|
|
|
351 |
chuckv |
388 |
|
352 |
mmeineke |
626 |
if (rijsq < rlistsq) then |
353 |
mmeineke |
377 |
|
354 |
|
|
nlist = nlist + 1 |
355 |
|
|
|
356 |
|
|
if (nlist > neighborListSize) then |
357 |
|
|
call expandNeighborList(natoms, listerror) |
358 |
|
|
if (listerror /= 0) then |
359 |
|
|
error = -1 |
360 |
|
|
write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded." |
361 |
|
|
return |
362 |
|
|
end if |
363 |
|
|
neighborListSize = size(list) |
364 |
|
|
endif |
365 |
|
|
|
366 |
|
|
list(nlist) = j |
367 |
|
|
|
368 |
mmeineke |
626 |
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
369 |
chuckv |
441 |
u_l, A, f, t, pot) |
370 |
mmeineke |
626 |
|
371 |
mmeineke |
377 |
endif |
372 |
|
|
enddo inner |
373 |
|
|
enddo |
374 |
|
|
|
375 |
|
|
point(natoms) = nlist + 1 |
376 |
|
|
|
377 |
|
|
else !! (update) |
378 |
|
|
|
379 |
|
|
! use the list to find the neighbors |
380 |
|
|
do i = 1, natoms-1 |
381 |
|
|
JBEG = POINT(i) |
382 |
|
|
JEND = POINT(i+1) - 1 |
383 |
|
|
! check thiat molecule i has neighbors |
384 |
|
|
if (jbeg .le. jend) then |
385 |
|
|
|
386 |
|
|
do jnab = jbeg, jend |
387 |
|
|
j = list(jnab) |
388 |
|
|
|
389 |
|
|
call get_interatomic_vector(q(:,i), q(:,j), d, rijsq) |
390 |
|
|
call do_pair(i, j, rijsq, d, do_pot, do_stress, & |
391 |
chuckv |
441 |
u_l, A, f, t, pot) |
392 |
mmeineke |
377 |
|
393 |
|
|
enddo |
394 |
|
|
endif |
395 |
|
|
enddo |
396 |
|
|
endif |
397 |
|
|
|
398 |
|
|
#endif |
399 |
|
|
|
400 |
|
|
! phew, done with main loop. |
401 |
|
|
|
402 |
|
|
#ifdef IS_MPI |
403 |
|
|
!!distribute forces |
404 |
chuckv |
438 |
|
405 |
|
|
f_temp = 0.0_dp |
406 |
|
|
call scatter(f_Row,f_temp,plan_row3d) |
407 |
|
|
do i = 1,nlocal |
408 |
|
|
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
409 |
|
|
end do |
410 |
|
|
|
411 |
|
|
f_temp = 0.0_dp |
412 |
mmeineke |
377 |
call scatter(f_Col,f_temp,plan_col3d) |
413 |
|
|
do i = 1,nlocal |
414 |
|
|
f(1:3,i) = f(1:3,i) + f_temp(1:3,i) |
415 |
|
|
end do |
416 |
|
|
|
417 |
|
|
if (FF_UsesDirectionalAtoms() .and. SimUsesDirectionalAtoms()) then |
418 |
chuckv |
438 |
t_temp = 0.0_dp |
419 |
|
|
call scatter(t_Row,t_temp,plan_row3d) |
420 |
|
|
do i = 1,nlocal |
421 |
|
|
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
422 |
|
|
end do |
423 |
|
|
t_temp = 0.0_dp |
424 |
mmeineke |
377 |
call scatter(t_Col,t_temp,plan_col3d) |
425 |
|
|
|
426 |
|
|
do i = 1,nlocal |
427 |
|
|
t(1:3,i) = t(1:3,i) + t_temp(1:3,i) |
428 |
|
|
end do |
429 |
|
|
endif |
430 |
|
|
|
431 |
|
|
if (do_pot) then |
432 |
|
|
! scatter/gather pot_row into the members of my column |
433 |
|
|
call scatter(pot_Row, pot_Temp, plan_row) |
434 |
chuckv |
439 |
|
435 |
mmeineke |
377 |
! scatter/gather pot_local into all other procs |
436 |
|
|
! add resultant to get total pot |
437 |
|
|
do i = 1, nlocal |
438 |
|
|
pot_local = pot_local + pot_Temp(i) |
439 |
|
|
enddo |
440 |
chuckv |
439 |
|
441 |
|
|
pot_Temp = 0.0_DP |
442 |
mmeineke |
377 |
|
443 |
|
|
call scatter(pot_Col, pot_Temp, plan_col) |
444 |
|
|
do i = 1, nlocal |
445 |
|
|
pot_local = pot_local + pot_Temp(i) |
446 |
|
|
enddo |
447 |
chuckv |
439 |
|
448 |
mmeineke |
377 |
endif |
449 |
|
|
#endif |
450 |
|
|
|
451 |
|
|
if (FF_RequiresPostpairCalc() .and. SimRequiresPostpairCalc()) then |
452 |
|
|
|
453 |
|
|
if (FF_uses_RF .and. SimUsesRF()) then |
454 |
|
|
|
455 |
|
|
#ifdef IS_MPI |
456 |
|
|
call scatter(rf_Row,rf,plan_row3d) |
457 |
|
|
call scatter(rf_Col,rf_Temp,plan_col3d) |
458 |
|
|
do i = 1,nlocal |
459 |
|
|
rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i) |
460 |
|
|
end do |
461 |
|
|
#endif |
462 |
|
|
|
463 |
|
|
do i = 1, getNlocal() |
464 |
|
|
|
465 |
|
|
rfpot = 0.0_DP |
466 |
|
|
#ifdef IS_MPI |
467 |
|
|
me_i = atid_row(i) |
468 |
|
|
#else |
469 |
|
|
me_i = atid(i) |
470 |
|
|
#endif |
471 |
|
|
call getElementProperty(atypes, me_i, "is_DP", is_DP_i) |
472 |
|
|
if ( is_DP_i ) then |
473 |
|
|
call getElementProperty(atypes, me_i, "dipole_moment", mu_i) |
474 |
|
|
!! The reaction field needs to include a self contribution |
475 |
|
|
!! to the field: |
476 |
|
|
call accumulate_self_rf(i, mu_i, u_l) |
477 |
|
|
!! Get the reaction field contribution to the |
478 |
|
|
!! potential and torques: |
479 |
|
|
call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot) |
480 |
|
|
#ifdef IS_MPI |
481 |
|
|
pot_local = pot_local + rfpot |
482 |
|
|
#else |
483 |
|
|
pot = pot + rfpot |
484 |
|
|
|
485 |
|
|
#endif |
486 |
|
|
endif |
487 |
|
|
enddo |
488 |
|
|
endif |
489 |
|
|
endif |
490 |
|
|
|
491 |
|
|
|
492 |
|
|
#ifdef IS_MPI |
493 |
|
|
|
494 |
|
|
if (do_pot) then |
495 |
chuckv |
441 |
pot = pot + pot_local |
496 |
mmeineke |
377 |
!! we assume the c code will do the allreduce to get the total potential |
497 |
|
|
!! we could do it right here if we needed to... |
498 |
|
|
endif |
499 |
|
|
|
500 |
|
|
if (do_stress) then |
501 |
gezelter |
490 |
call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, & |
502 |
mmeineke |
377 |
mpi_comm_world,mpi_err) |
503 |
chuckv |
470 |
call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, & |
504 |
mmeineke |
377 |
mpi_comm_world,mpi_err) |
505 |
|
|
endif |
506 |
|
|
|
507 |
|
|
#else |
508 |
|
|
|
509 |
|
|
if (do_stress) then |
510 |
|
|
tau = tau_Temp |
511 |
|
|
virial = virial_Temp |
512 |
|
|
endif |
513 |
|
|
|
514 |
|
|
#endif |
515 |
|
|
|
516 |
|
|
end subroutine do_force_loop |
517 |
|
|
|
518 |
chuckv |
441 |
subroutine do_pair(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot) |
519 |
mmeineke |
377 |
|
520 |
|
|
real( kind = dp ) :: pot |
521 |
chuckv |
460 |
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
522 |
|
|
real (kind=dp), dimension(9,getNlocal()) :: A |
523 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: f |
524 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: t |
525 |
mmeineke |
377 |
|
526 |
|
|
logical, intent(inout) :: do_pot, do_stress |
527 |
|
|
integer, intent(in) :: i, j |
528 |
|
|
real ( kind = dp ), intent(inout) :: rijsq |
529 |
|
|
real ( kind = dp ) :: r |
530 |
|
|
real ( kind = dp ), intent(inout) :: d(3) |
531 |
|
|
logical :: is_LJ_i, is_LJ_j |
532 |
|
|
logical :: is_DP_i, is_DP_j |
533 |
|
|
logical :: is_GB_i, is_GB_j |
534 |
|
|
logical :: is_Sticky_i, is_Sticky_j |
535 |
|
|
integer :: me_i, me_j |
536 |
|
|
|
537 |
|
|
r = sqrt(rijsq) |
538 |
|
|
|
539 |
|
|
#ifdef IS_MPI |
540 |
gezelter |
490 |
if (tagRow(i) .eq. tagColumn(j)) then |
541 |
|
|
write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j) |
542 |
|
|
endif |
543 |
mmeineke |
377 |
|
544 |
|
|
me_i = atid_row(i) |
545 |
|
|
me_j = atid_col(j) |
546 |
|
|
|
547 |
|
|
#else |
548 |
|
|
|
549 |
|
|
me_i = atid(i) |
550 |
|
|
me_j = atid(j) |
551 |
|
|
|
552 |
|
|
#endif |
553 |
|
|
|
554 |
|
|
if (FF_uses_LJ .and. SimUsesLJ()) then |
555 |
|
|
call getElementProperty(atypes, me_i, "is_LJ", is_LJ_i) |
556 |
|
|
call getElementProperty(atypes, me_j, "is_LJ", is_LJ_j) |
557 |
|
|
|
558 |
|
|
if ( is_LJ_i .and. is_LJ_j ) & |
559 |
|
|
call do_lj_pair(i, j, d, r, rijsq, pot, f, do_pot, do_stress) |
560 |
|
|
endif |
561 |
|
|
|
562 |
|
|
if (FF_uses_dipoles .and. SimUsesDipoles()) then |
563 |
|
|
call getElementProperty(atypes, me_i, "is_DP", is_DP_i) |
564 |
|
|
call getElementProperty(atypes, me_j, "is_DP", is_DP_j) |
565 |
|
|
|
566 |
|
|
if ( is_DP_i .and. is_DP_j ) then |
567 |
|
|
|
568 |
gezelter |
462 |
call do_dipole_pair(i, j, d, r, rijsq, pot, u_l, f, t, & |
569 |
mmeineke |
377 |
do_pot, do_stress) |
570 |
|
|
if (FF_uses_RF .and. SimUsesRF()) then |
571 |
|
|
call accumulate_rf(i, j, r, u_l) |
572 |
|
|
call rf_correct_forces(i, j, d, r, u_l, f, do_stress) |
573 |
|
|
endif |
574 |
|
|
|
575 |
|
|
endif |
576 |
|
|
endif |
577 |
|
|
|
578 |
|
|
if (FF_uses_Sticky .and. SimUsesSticky()) then |
579 |
|
|
|
580 |
|
|
call getElementProperty(atypes, me_i, "is_Sticky", is_Sticky_i) |
581 |
|
|
call getElementProperty(atypes, me_j, "is_Sticky", is_Sticky_j) |
582 |
chuckv |
388 |
|
583 |
mmeineke |
377 |
if ( is_Sticky_i .and. is_Sticky_j ) then |
584 |
|
|
call do_sticky_pair(i, j, d, r, rijsq, A, pot, f, t, & |
585 |
|
|
do_pot, do_stress) |
586 |
|
|
endif |
587 |
|
|
endif |
588 |
|
|
|
589 |
|
|
|
590 |
|
|
if (FF_uses_GB .and. SimUsesGB()) then |
591 |
|
|
|
592 |
|
|
call getElementProperty(atypes, me_i, "is_GB", is_GB_i) |
593 |
|
|
call getElementProperty(atypes, me_j, "is_GB", is_GB_j) |
594 |
|
|
|
595 |
|
|
if ( is_GB_i .and. is_GB_j ) then |
596 |
|
|
call do_gb_pair(i, j, d, r, rijsq, u_l, pot, f, t, & |
597 |
|
|
do_pot, do_stress) |
598 |
|
|
endif |
599 |
|
|
endif |
600 |
|
|
|
601 |
mmeineke |
597 |
|
602 |
|
|
|
603 |
mmeineke |
377 |
end subroutine do_pair |
604 |
|
|
|
605 |
|
|
|
606 |
chuckv |
631 |
|
607 |
|
|
subroutine do_preforce(i, j, rijsq, d, do_pot, do_stress, u_l, A, f, t, pot) |
608 |
|
|
real( kind = dp ) :: pot |
609 |
|
|
real( kind = dp ), dimension(3,getNlocal()) :: u_l |
610 |
|
|
real (kind=dp), dimension(9,getNlocal()) :: A |
611 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: f |
612 |
|
|
real (kind=dp), dimension(3,getNlocal()) :: t |
613 |
|
|
|
614 |
|
|
logical, intent(inout) :: do_pot, do_stress |
615 |
|
|
integer, intent(in) :: i, j |
616 |
|
|
real ( kind = dp ), intent(inout) :: rijsq |
617 |
|
|
real ( kind = dp ) :: r |
618 |
|
|
real ( kind = dp ), intent(inout) :: d(3) |
619 |
|
|
|
620 |
|
|
logical :: is_EAM_i, is_EAM_j |
621 |
|
|
|
622 |
|
|
integer :: me_i, me_j |
623 |
|
|
|
624 |
|
|
r = sqrt(rijsq) |
625 |
|
|
|
626 |
|
|
#ifdef IS_MPI |
627 |
|
|
if (tagRow(i) .eq. tagColumn(j)) then |
628 |
|
|
write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j) |
629 |
|
|
endif |
630 |
|
|
|
631 |
|
|
me_i = atid_row(i) |
632 |
|
|
me_j = atid_col(j) |
633 |
|
|
|
634 |
|
|
#else |
635 |
|
|
|
636 |
|
|
me_i = atid(i) |
637 |
|
|
me_j = atid(j) |
638 |
|
|
|
639 |
|
|
#endif |
640 |
|
|
|
641 |
|
|
if (FF_uses_EAM .and. SimUsesEAM()) then |
642 |
|
|
call getElementProperty(atypes, me_i, "is_EAM", is_EAM_i) |
643 |
|
|
call getElementProperty(atypes, me_j, "is_EAM", is_EAM_j) |
644 |
|
|
|
645 |
mmeineke |
634 |
!!$ if ( is_EAM_i .and. is_EAM_j ) & |
646 |
|
|
!!$ call calc_EAM_prepair(i, j, d, r, rijsq ) |
647 |
chuckv |
631 |
endif |
648 |
|
|
end subroutine do_preforce |
649 |
|
|
|
650 |
|
|
|
651 |
mmeineke |
377 |
subroutine get_interatomic_vector(q_i, q_j, d, r_sq) |
652 |
|
|
|
653 |
|
|
real (kind = dp), dimension(3) :: q_i |
654 |
|
|
real (kind = dp), dimension(3) :: q_j |
655 |
|
|
real ( kind = dp ), intent(out) :: r_sq |
656 |
gezelter |
571 |
real( kind = dp ) :: d(3), scaled(3) |
657 |
|
|
integer i |
658 |
|
|
|
659 |
chuckv |
482 |
d(1:3) = q_j(1:3) - q_i(1:3) |
660 |
gezelter |
571 |
|
661 |
mmeineke |
377 |
! Wrap back into periodic box if necessary |
662 |
|
|
if ( SimUsesPBC() ) then |
663 |
mmeineke |
393 |
|
664 |
gezelter |
571 |
if( .not.boxIsOrthorhombic ) then |
665 |
|
|
! calc the scaled coordinates. |
666 |
|
|
|
667 |
mmeineke |
572 |
scaled = matmul(HmatInv, d) |
668 |
gezelter |
571 |
|
669 |
|
|
! wrap the scaled coordinates |
670 |
|
|
|
671 |
mmeineke |
572 |
scaled = scaled - anint(scaled) |
672 |
|
|
|
673 |
gezelter |
571 |
|
674 |
|
|
! calc the wrapped real coordinates from the wrapped scaled |
675 |
|
|
! coordinates |
676 |
|
|
|
677 |
mmeineke |
572 |
d = matmul(Hmat,scaled) |
678 |
gezelter |
571 |
|
679 |
|
|
else |
680 |
|
|
! calc the scaled coordinates. |
681 |
|
|
|
682 |
|
|
do i = 1, 3 |
683 |
|
|
scaled(i) = d(i) * HmatInv(i,i) |
684 |
|
|
|
685 |
|
|
! wrap the scaled coordinates |
686 |
|
|
|
687 |
|
|
scaled(i) = scaled(i) - anint(scaled(i)) |
688 |
|
|
|
689 |
|
|
! calc the wrapped real coordinates from the wrapped scaled |
690 |
|
|
! coordinates |
691 |
|
|
|
692 |
|
|
d(i) = scaled(i)*Hmat(i,i) |
693 |
|
|
enddo |
694 |
|
|
endif |
695 |
mmeineke |
393 |
|
696 |
mmeineke |
377 |
endif |
697 |
gezelter |
571 |
|
698 |
mmeineke |
377 |
r_sq = dot_product(d,d) |
699 |
gezelter |
571 |
|
700 |
mmeineke |
377 |
end subroutine get_interatomic_vector |
701 |
gezelter |
571 |
|
702 |
mmeineke |
377 |
subroutine check_initialization(error) |
703 |
|
|
integer, intent(out) :: error |
704 |
|
|
|
705 |
|
|
error = 0 |
706 |
|
|
! Make sure we are properly initialized. |
707 |
|
|
if (.not. do_forces_initialized) then |
708 |
|
|
error = -1 |
709 |
|
|
return |
710 |
|
|
endif |
711 |
|
|
|
712 |
|
|
#ifdef IS_MPI |
713 |
|
|
if (.not. isMPISimSet()) then |
714 |
|
|
write(default_error,*) "ERROR: mpiSimulation has not been initialized!" |
715 |
|
|
error = -1 |
716 |
|
|
return |
717 |
|
|
endif |
718 |
|
|
#endif |
719 |
|
|
|
720 |
|
|
return |
721 |
|
|
end subroutine check_initialization |
722 |
|
|
|
723 |
|
|
|
724 |
|
|
subroutine zero_work_arrays() |
725 |
|
|
|
726 |
|
|
#ifdef IS_MPI |
727 |
|
|
|
728 |
|
|
q_Row = 0.0_dp |
729 |
|
|
q_Col = 0.0_dp |
730 |
|
|
|
731 |
|
|
u_l_Row = 0.0_dp |
732 |
|
|
u_l_Col = 0.0_dp |
733 |
|
|
|
734 |
|
|
A_Row = 0.0_dp |
735 |
|
|
A_Col = 0.0_dp |
736 |
|
|
|
737 |
|
|
f_Row = 0.0_dp |
738 |
|
|
f_Col = 0.0_dp |
739 |
|
|
f_Temp = 0.0_dp |
740 |
|
|
|
741 |
|
|
t_Row = 0.0_dp |
742 |
|
|
t_Col = 0.0_dp |
743 |
|
|
t_Temp = 0.0_dp |
744 |
|
|
|
745 |
|
|
pot_Row = 0.0_dp |
746 |
|
|
pot_Col = 0.0_dp |
747 |
|
|
pot_Temp = 0.0_dp |
748 |
|
|
|
749 |
|
|
rf_Row = 0.0_dp |
750 |
|
|
rf_Col = 0.0_dp |
751 |
|
|
rf_Temp = 0.0_dp |
752 |
|
|
|
753 |
|
|
#endif |
754 |
|
|
|
755 |
|
|
rf = 0.0_dp |
756 |
|
|
tau_Temp = 0.0_dp |
757 |
|
|
virial_Temp = 0.0_dp |
758 |
|
|
end subroutine zero_work_arrays |
759 |
|
|
|
760 |
|
|
function skipThisPair(atom1, atom2) result(skip_it) |
761 |
|
|
integer, intent(in) :: atom1 |
762 |
|
|
integer, intent(in), optional :: atom2 |
763 |
|
|
logical :: skip_it |
764 |
|
|
integer :: unique_id_1, unique_id_2 |
765 |
chuckv |
388 |
integer :: me_i,me_j |
766 |
mmeineke |
377 |
integer :: i |
767 |
|
|
|
768 |
|
|
skip_it = .false. |
769 |
|
|
|
770 |
|
|
!! there are a number of reasons to skip a pair or a particle |
771 |
|
|
!! mostly we do this to exclude atoms who are involved in short |
772 |
|
|
!! range interactions (bonds, bends, torsions), but we also need |
773 |
|
|
!! to exclude some overcounted interactions that result from |
774 |
|
|
!! the parallel decomposition |
775 |
|
|
|
776 |
|
|
#ifdef IS_MPI |
777 |
|
|
!! in MPI, we have to look up the unique IDs for each atom |
778 |
|
|
unique_id_1 = tagRow(atom1) |
779 |
|
|
#else |
780 |
|
|
!! in the normal loop, the atom numbers are unique |
781 |
|
|
unique_id_1 = atom1 |
782 |
|
|
#endif |
783 |
chuckv |
388 |
|
784 |
mmeineke |
377 |
!! We were called with only one atom, so just check the global exclude |
785 |
|
|
!! list for this atom |
786 |
|
|
if (.not. present(atom2)) then |
787 |
|
|
do i = 1, nExcludes_global |
788 |
|
|
if (excludesGlobal(i) == unique_id_1) then |
789 |
|
|
skip_it = .true. |
790 |
|
|
return |
791 |
|
|
end if |
792 |
|
|
end do |
793 |
|
|
return |
794 |
|
|
end if |
795 |
|
|
|
796 |
|
|
#ifdef IS_MPI |
797 |
|
|
unique_id_2 = tagColumn(atom2) |
798 |
|
|
#else |
799 |
|
|
unique_id_2 = atom2 |
800 |
|
|
#endif |
801 |
chuckv |
441 |
|
802 |
mmeineke |
377 |
#ifdef IS_MPI |
803 |
|
|
!! this situation should only arise in MPI simulations |
804 |
|
|
if (unique_id_1 == unique_id_2) then |
805 |
|
|
skip_it = .true. |
806 |
|
|
return |
807 |
|
|
end if |
808 |
|
|
|
809 |
|
|
!! this prevents us from doing the pair on multiple processors |
810 |
|
|
if (unique_id_1 < unique_id_2) then |
811 |
chuckv |
441 |
if (mod(unique_id_1 + unique_id_2,2) == 0) then |
812 |
|
|
skip_it = .true. |
813 |
|
|
return |
814 |
|
|
endif |
815 |
mmeineke |
377 |
else |
816 |
chuckv |
441 |
if (mod(unique_id_1 + unique_id_2,2) == 1) then |
817 |
|
|
skip_it = .true. |
818 |
|
|
return |
819 |
|
|
endif |
820 |
mmeineke |
377 |
endif |
821 |
|
|
#endif |
822 |
chuckv |
441 |
|
823 |
mmeineke |
377 |
!! the rest of these situations can happen in all simulations: |
824 |
|
|
do i = 1, nExcludes_global |
825 |
|
|
if ((excludesGlobal(i) == unique_id_1) .or. & |
826 |
|
|
(excludesGlobal(i) == unique_id_2)) then |
827 |
|
|
skip_it = .true. |
828 |
|
|
return |
829 |
|
|
endif |
830 |
|
|
enddo |
831 |
chuckv |
441 |
|
832 |
mmeineke |
377 |
do i = 1, nExcludes_local |
833 |
|
|
if (excludesLocal(1,i) == unique_id_1) then |
834 |
|
|
if (excludesLocal(2,i) == unique_id_2) then |
835 |
|
|
skip_it = .true. |
836 |
|
|
return |
837 |
|
|
endif |
838 |
|
|
else |
839 |
|
|
if (excludesLocal(1,i) == unique_id_2) then |
840 |
|
|
if (excludesLocal(2,i) == unique_id_1) then |
841 |
|
|
skip_it = .true. |
842 |
|
|
return |
843 |
|
|
endif |
844 |
|
|
endif |
845 |
|
|
endif |
846 |
|
|
end do |
847 |
|
|
|
848 |
|
|
return |
849 |
|
|
end function skipThisPair |
850 |
|
|
|
851 |
|
|
function FF_UsesDirectionalAtoms() result(doesit) |
852 |
|
|
logical :: doesit |
853 |
|
|
doesit = FF_uses_dipoles .or. FF_uses_sticky .or. & |
854 |
|
|
FF_uses_GB .or. FF_uses_RF |
855 |
|
|
end function FF_UsesDirectionalAtoms |
856 |
|
|
|
857 |
|
|
function FF_RequiresPrepairCalc() result(doesit) |
858 |
|
|
logical :: doesit |
859 |
|
|
doesit = FF_uses_EAM |
860 |
|
|
end function FF_RequiresPrepairCalc |
861 |
|
|
|
862 |
|
|
function FF_RequiresPostpairCalc() result(doesit) |
863 |
|
|
logical :: doesit |
864 |
|
|
doesit = FF_uses_RF |
865 |
|
|
end function FF_RequiresPostpairCalc |
866 |
|
|
|
867 |
|
|
end module do_Forces |