ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 597
Committed: Mon Jul 14 21:28:54 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 22811 byte(s)
Log Message:
found a bug. Unit vectors were not being updated

File Contents

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