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