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