ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 1169
Committed: Wed May 12 19:44:54 2004 UTC (20 years, 4 months ago) by gezelter
File size: 45481 byte(s)
Log Message:
MPI fixes and removal of extraneous write statements

File Contents

# User Rev Content
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 gezelter 1169 !! @version $Id: do_Forces.F90,v 1.58 2004-05-12 19:44:49 gezelter Exp $, $Date: 2004-05-12 19:44:49 $, $Name: not supported by cvs2svn $, $Revision: 1.58 $
8 mmeineke 377
9     module do_Forces
10     use force_globals
11     use simulation
12     use definitions
13     use atype_module
14 gezelter 1150 use switcheroo
15 mmeineke 377 use neighborLists
16     use lj
17     use sticky_pair
18     use dipole_dipole
19 gezelter 941 use charge_charge
20 mmeineke 377 use reaction_field
21     use gb_pair
22 mmeineke 626 use vector_class
23 chuckv 650 use eam
24 chuckv 657 use status
25 mmeineke 377 #ifdef IS_MPI
26     use mpiSimulation
27     #endif
28    
29     implicit none
30     PRIVATE
31    
32     #define __FORTRAN90
33     #include "fForceField.h"
34 gezelter 1150 #include "fSwitchingFunction.h"
35 mmeineke 377
36 chuckv 900 logical, save :: haveRlist = .false.
37     logical, save :: haveNeighborList = .false.
38 mmeineke 626 logical, save :: havePolicies = .false.
39 chuckv 900 logical, save :: haveSIMvariables = .false.
40     logical, save :: havePropertyMap = .false.
41     logical, save :: haveSaneForceField = .false.
42 mmeineke 377 logical, save :: FF_uses_LJ
43     logical, save :: FF_uses_sticky
44 gezelter 941 logical, save :: FF_uses_charges
45 mmeineke 377 logical, save :: FF_uses_dipoles
46     logical, save :: FF_uses_RF
47     logical, save :: FF_uses_GB
48     logical, save :: FF_uses_EAM
49 chuckv 900 logical, save :: SIM_uses_LJ
50     logical, save :: SIM_uses_sticky
51 gezelter 941 logical, save :: SIM_uses_charges
52 chuckv 900 logical, save :: SIM_uses_dipoles
53     logical, save :: SIM_uses_RF
54     logical, save :: SIM_uses_GB
55     logical, save :: SIM_uses_EAM
56     logical, save :: SIM_requires_postpair_calc
57     logical, save :: SIM_requires_prepair_calc
58     logical, save :: SIM_uses_directional_atoms
59     logical, save :: SIM_uses_PBC
60 gezelter 1138 logical, save :: SIM_uses_molecular_cutoffs
61 mmeineke 377
62 mmeineke 626 real(kind=dp), save :: rlist, rlistsq
63    
64 mmeineke 377 public :: init_FF
65     public :: do_force_loop
66 mmeineke 626 public :: setRlistDF
67 mmeineke 377
68 chuckv 694 #ifdef PROFILE
69 chuckv 883 public :: getforcetime
70     real, save :: forceTime = 0
71     real :: forceTimeInitial, forceTimeFinal
72 mmeineke 891 integer :: nLoops
73 chuckv 694 #endif
74    
75 chuckv 900 type :: Properties
76     logical :: is_lj = .false.
77     logical :: is_sticky = .false.
78     logical :: is_dp = .false.
79     logical :: is_gb = .false.
80     logical :: is_eam = .false.
81 gezelter 941 logical :: is_charge = .false.
82     real(kind=DP) :: charge = 0.0_DP
83 chuckv 900 real(kind=DP) :: dipole_moment = 0.0_DP
84     end type Properties
85 chuckv 895
86 chuckv 900 type(Properties), dimension(:),allocatable :: PropertyMap
87    
88 mmeineke 377 contains
89    
90 mmeineke 626 subroutine setRlistDF( this_rlist )
91    
92     real(kind=dp) :: this_rlist
93    
94     rlist = this_rlist
95     rlistsq = rlist * rlist
96    
97     haveRlist = .true.
98    
99     end subroutine setRlistDF
100    
101 chuckv 900 subroutine createPropertyMap(status)
102     integer :: nAtypes
103     integer :: status
104     integer :: i
105     logical :: thisProperty
106     real (kind=DP) :: thisDPproperty
107    
108     status = 0
109    
110     nAtypes = getSize(atypes)
111    
112     if (nAtypes == 0) then
113     status = -1
114     return
115     end if
116    
117     if (.not. allocated(PropertyMap)) then
118     allocate(PropertyMap(nAtypes))
119     endif
120    
121     do i = 1, nAtypes
122     call getElementProperty(atypes, i, "is_LJ", thisProperty)
123     PropertyMap(i)%is_LJ = thisProperty
124 gezelter 941
125     call getElementProperty(atypes, i, "is_Charge", thisProperty)
126     PropertyMap(i)%is_Charge = thisProperty
127    
128     if (thisProperty) then
129     call getElementProperty(atypes, i, "charge", thisDPproperty)
130     PropertyMap(i)%charge = thisDPproperty
131     endif
132    
133 chuckv 900 call getElementProperty(atypes, i, "is_DP", thisProperty)
134     PropertyMap(i)%is_DP = thisProperty
135    
136     if (thisProperty) then
137     call getElementProperty(atypes, i, "dipole_moment", thisDPproperty)
138     PropertyMap(i)%dipole_moment = thisDPproperty
139     endif
140    
141     call getElementProperty(atypes, i, "is_Sticky", thisProperty)
142     PropertyMap(i)%is_Sticky = thisProperty
143     call getElementProperty(atypes, i, "is_GB", thisProperty)
144     PropertyMap(i)%is_GB = thisProperty
145     call getElementProperty(atypes, i, "is_EAM", thisProperty)
146     PropertyMap(i)%is_EAM = thisProperty
147     end do
148    
149     havePropertyMap = .true.
150    
151     end subroutine createPropertyMap
152    
153     subroutine setSimVariables()
154     SIM_uses_LJ = SimUsesLJ()
155     SIM_uses_sticky = SimUsesSticky()
156 gezelter 941 SIM_uses_charges = SimUsesCharges()
157 chuckv 900 SIM_uses_dipoles = SimUsesDipoles()
158     SIM_uses_RF = SimUsesRF()
159     SIM_uses_GB = SimUsesGB()
160     SIM_uses_EAM = SimUsesEAM()
161     SIM_requires_postpair_calc = SimRequiresPostpairCalc()
162     SIM_requires_prepair_calc = SimRequiresPrepairCalc()
163     SIM_uses_directional_atoms = SimUsesDirectionalAtoms()
164     SIM_uses_PBC = SimUsesPBC()
165 tim 1144 !SIM_uses_molecular_cutoffs = SimUsesMolecularCutoffs()
166 chuckv 900
167     haveSIMvariables = .true.
168    
169     return
170     end subroutine setSimVariables
171    
172     subroutine doReadyCheck(error)
173     integer, intent(out) :: error
174    
175     integer :: myStatus
176    
177     error = 0
178    
179     if (.not. havePropertyMap) then
180    
181     myStatus = 0
182    
183     call createPropertyMap(myStatus)
184    
185     if (myStatus .ne. 0) then
186     write(default_error, *) 'createPropertyMap failed in do_Forces!'
187     error = -1
188     return
189     endif
190     endif
191    
192     if (.not. haveSIMvariables) then
193     call setSimVariables()
194     endif
195    
196     if (.not. haveRlist) then
197     write(default_error, *) 'rList has not been set in do_Forces!'
198     error = -1
199     return
200     endif
201    
202     if (SIM_uses_LJ .and. FF_uses_LJ) then
203     if (.not. havePolicies) then
204     write(default_error, *) 'LJ mixing Policies have not been set in do_Forces!'
205     error = -1
206     return
207     endif
208     endif
209    
210     if (.not. haveNeighborList) then
211     write(default_error, *) 'neighbor list has not been initialized in do_Forces!'
212     error = -1
213     return
214     end if
215    
216     if (.not. haveSaneForceField) then
217     write(default_error, *) 'Force Field is not sane in do_Forces!'
218     error = -1
219     return
220     end if
221    
222     #ifdef IS_MPI
223     if (.not. isMPISimSet()) then
224     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
225     error = -1
226     return
227     endif
228     #endif
229     return
230     end subroutine doReadyCheck
231    
232    
233 mmeineke 377 subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat)
234    
235     integer, intent(in) :: LJMIXPOLICY
236     logical, intent(in) :: use_RF_c
237    
238     integer, intent(out) :: thisStat
239     integer :: my_status, nMatches
240     integer, pointer :: MatchList(:) => null()
241     real(kind=dp) :: rcut, rrf, rt, dielect
242    
243     !! assume things are copacetic, unless they aren't
244     thisStat = 0
245    
246     !! Fortran's version of a cast:
247     FF_uses_RF = use_RF_c
248    
249     !! init_FF is called *after* all of the atom types have been
250     !! defined in atype_module using the new_atype subroutine.
251     !!
252     !! this will scan through the known atypes and figure out what
253     !! interactions are used by the force field.
254    
255     FF_uses_LJ = .false.
256     FF_uses_sticky = .false.
257 gezelter 941 FF_uses_charges = .false.
258 mmeineke 377 FF_uses_dipoles = .false.
259     FF_uses_GB = .false.
260     FF_uses_EAM = .false.
261    
262     call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList)
263     if (nMatches .gt. 0) FF_uses_LJ = .true.
264 gezelter 941
265     call getMatchingElementList(atypes, "is_Charge", .true., nMatches, MatchList)
266 tim 1113 if (nMatches .gt. 0) FF_uses_charges = .true.
267    
268 mmeineke 377 call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList)
269     if (nMatches .gt. 0) FF_uses_dipoles = .true.
270    
271     call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, &
272     MatchList)
273     if (nMatches .gt. 0) FF_uses_Sticky = .true.
274    
275     call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList)
276     if (nMatches .gt. 0) FF_uses_GB = .true.
277    
278     call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
279     if (nMatches .gt. 0) FF_uses_EAM = .true.
280    
281 chuckv 900 !! Assume sanity (for the sake of argument)
282     haveSaneForceField = .true.
283    
284 mmeineke 377 !! check to make sure the FF_uses_RF setting makes sense
285    
286     if (FF_uses_dipoles) then
287     if (FF_uses_RF) then
288     dielect = getDielect()
289 mmeineke 626 call initialize_rf(dielect)
290 mmeineke 377 endif
291     else
292     if (FF_uses_RF) then
293     write(default_error,*) 'Using Reaction Field with no dipoles? Huh?'
294     thisStat = -1
295 chuckv 900 haveSaneForceField = .false.
296 mmeineke 377 return
297     endif
298 mmeineke 626 endif
299 mmeineke 377
300     if (FF_uses_LJ) then
301    
302     select case (LJMIXPOLICY)
303 gezelter 834 case (LB_MIXING_RULE)
304     call init_lj_FF(LB_MIXING_RULE, my_status)
305     case (EXPLICIT_MIXING_RULE)
306     call init_lj_FF(EXPLICIT_MIXING_RULE, my_status)
307 mmeineke 377 case default
308     write(default_error,*) 'unknown LJ Mixing Policy!'
309     thisStat = -1
310 chuckv 900 haveSaneForceField = .false.
311 mmeineke 377 return
312     end select
313     if (my_status /= 0) then
314     thisStat = -1
315 chuckv 900 haveSaneForceField = .false.
316 mmeineke 377 return
317     end if
318 chuckv 900 havePolicies = .true.
319 mmeineke 377 endif
320    
321     if (FF_uses_sticky) then
322     call check_sticky_FF(my_status)
323     if (my_status /= 0) then
324     thisStat = -1
325 chuckv 900 haveSaneForceField = .false.
326 mmeineke 377 return
327     end if
328     endif
329 chuckv 657
330    
331     if (FF_uses_EAM) then
332 chuckv 801 call init_EAM_FF(my_status)
333 chuckv 657 if (my_status /= 0) then
334 chuckv 900 write(default_error, *) "init_EAM_FF returned a bad status"
335 chuckv 657 thisStat = -1
336 chuckv 900 haveSaneForceField = .false.
337 chuckv 657 return
338     end if
339     endif
340    
341 mmeineke 377 if (FF_uses_GB) then
342     call check_gb_pair_FF(my_status)
343     if (my_status .ne. 0) then
344     thisStat = -1
345 chuckv 900 haveSaneForceField = .false.
346 mmeineke 377 return
347     endif
348     endif
349    
350     if (FF_uses_GB .and. FF_uses_LJ) then
351     endif
352 chuckv 900 if (.not. haveNeighborList) then
353 chuckv 480 !! Create neighbor lists
354 chuckv 898 call expandNeighborList(nLocal, my_status)
355 chuckv 480 if (my_Status /= 0) then
356     write(default_error,*) "SimSetup: ExpandNeighborList returned error."
357     thisStat = -1
358     return
359     endif
360 chuckv 900 haveNeighborList = .true.
361 chuckv 480 endif
362 gezelter 1150
363 chuckv 657
364 gezelter 1150
365 mmeineke 377 end subroutine init_FF
366    
367    
368     !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
369     !------------------------------------------------------------->
370 gezelter 1150 subroutine do_force_loop(q, q_group, A, u_l, f, t, tau, pot, &
371 gezelter 1138 do_pot_c, do_stress_c, error)
372 mmeineke 377 !! Position array provided by C, dimensioned by getNlocal
373 gezelter 1166 real ( kind = dp ), dimension(3, nLocal) :: q
374 gezelter 1138 !! molecular center-of-mass position array
375 gezelter 1166 real ( kind = dp ), dimension(3, nGroup) :: q_group
376 mmeineke 377 !! Rotation Matrix for each long range particle in simulation.
377 gezelter 1166 real( kind = dp), dimension(9, nLocal) :: A
378 mmeineke 377 !! Unit vectors for dipoles (lab frame)
379 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
380 mmeineke 377 !! Force array provided by C, dimensioned by getNlocal
381 chuckv 898 real ( kind = dp ), dimension(3,nLocal) :: f
382 mmeineke 377 !! Torsion array provided by C, dimensioned by getNlocal
383 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: t
384 chuckv 895
385 mmeineke 377 !! Stress Tensor
386     real( kind = dp), dimension(9) :: tau
387     real ( kind = dp ) :: pot
388     logical ( kind = 2) :: do_pot_c, do_stress_c
389     logical :: do_pot
390     logical :: do_stress
391 gezelter 1150 logical :: in_switching_region
392 chuckv 439 #ifdef IS_MPI
393 chuckv 441 real( kind = DP ) :: pot_local
394 mmeineke 377 integer :: nrow
395     integer :: ncol
396 chuckv 694 integer :: nprocs
397 gezelter 1150 integer :: nrow_group
398     integer :: ncol_group
399 mmeineke 377 #endif
400     integer :: natoms
401     logical :: update_nlist
402     integer :: i, j, jbeg, jend, jnab
403 gezelter 1150 integer :: ia, jb, atom1, atom2
404 mmeineke 377 integer :: nlist
405 gezelter 1169 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, vpair, vij
406 gezelter 1150 real( kind = DP ) :: sw, dswdr, swderiv, mf
407     real(kind=dp),dimension(3) :: d_atm, d_grp
408 mmeineke 377 real(kind=dp) :: rfpot, mu_i, virial
409 gezelter 1169 integer :: me_i, me_j, n_in_i, n_in_j
410 mmeineke 377 logical :: is_dp_i
411     integer :: neighborListSize
412     integer :: listerror, error
413     integer :: localError
414 chuckv 897 integer :: propPack_i, propPack_j
415 mmeineke 626
416 gezelter 845 real(kind=dp) :: listSkin = 1.0
417 gezelter 1138
418 mmeineke 377 !! initialize local variables
419 gezelter 1138
420 mmeineke 377 #ifdef IS_MPI
421 chuckv 441 pot_local = 0.0_dp
422 mmeineke 377 nrow = getNrow(plan_row)
423     ncol = getNcol(plan_col)
424 gezelter 1150 nrow_group = getNrowGroup(plan_row)
425     ncol_group = getNcolGroup(plan_col)
426 mmeineke 377 #else
427     natoms = nlocal
428     #endif
429 gezelter 1138
430 chuckv 900 call doReadyCheck(localError)
431 mmeineke 377 if ( localError .ne. 0 ) then
432 chuckv 900 call handleError("do_force_loop", "Not Initialized")
433 mmeineke 377 error = -1
434     return
435     end if
436     call zero_work_arrays()
437 gezelter 1138
438 mmeineke 377 do_pot = do_pot_c
439     do_stress = do_stress_c
440 gezelter 1138
441 mmeineke 377 ! Gather all information needed by all force loops:
442    
443     #ifdef IS_MPI
444 gezelter 1138
445 gezelter 1150 call gather(q, q_Row, plan_row3d)
446     call gather(q, q_Col, plan_col3d)
447    
448     call gather(q_group, q_group_Row, plan_row_Group_3d)
449     call gather(q_group, q_group_Col, plan_col_Group_3d)
450    
451 chuckv 900 if (FF_UsesDirectionalAtoms() .and. SIM_uses_directional_atoms) then
452 mmeineke 377 call gather(u_l,u_l_Row,plan_row3d)
453     call gather(u_l,u_l_Col,plan_col3d)
454    
455     call gather(A,A_Row,plan_row_rotation)
456     call gather(A,A_Col,plan_col_rotation)
457     endif
458    
459     #endif
460 gezelter 1138
461     !! Begin force loop timing:
462 chuckv 694 #ifdef PROFILE
463     call cpu_time(forceTimeInitial)
464     nloops = nloops + 1
465     #endif
466 gezelter 1138
467 chuckv 900 if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
468 mmeineke 377 !! See if we need to update neighbor lists
469 gezelter 1150
470     call checkNeighborList(nGroup, q_group, listSkin, update_nlist)
471    
472 mmeineke 377 !! if_mpi_gather_stuff_for_prepair
473     !! do_prepair_loop_if_needed
474     !! if_mpi_scatter_stuff_from_prepair
475     !! if_mpi_gather_stuff_from_prepair_to_main_loop
476 gezelter 1138
477     !--------------------PREFORCE LOOP----------->>>>>>>>>>>>>>>>>>>>>>>>>>>
478 chuckv 648 #ifdef IS_MPI
479    
480 gezelter 1138 if (update_nlist) then
481 chuckv 648
482 gezelter 1138 !! save current configuration, construct neighbor list,
483     !! and calculate forces
484 gezelter 1150
485     call saveNeighborList(nGroup, q_group)
486 gezelter 1138
487     neighborListSize = size(list)
488     nlist = 0
489    
490 gezelter 1150 do i = 1, nrow_group
491 gezelter 1138 point(i) = nlist + 1
492 chuckv 648
493 gezelter 1150 do j = 1, ncol_group
494 chuckv 648
495 gezelter 1150 call get_interatomic_vector(q_group_Row(:,i), &
496     q_group_Col(:,j), d_grp, rgrpsq)
497 chuckv 648
498 gezelter 1150 if (rgrpsq < rlistsq) then
499 gezelter 1138 nlist = nlist + 1
500    
501     if (nlist > neighborListSize) then
502 gezelter 1150 call expandNeighborList(nGroup, listerror)
503 gezelter 1138 if (listerror /= 0) then
504     error = -1
505     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
506     return
507     end if
508     neighborListSize = size(list)
509     endif
510    
511     list(nlist) = j
512    
513 gezelter 1150 do ia = groupStart(i), groupStart(i+1)-1
514     atom1 = groupList(ia)
515    
516     prepair_inner: do jb = groupStart(j), groupStart(j+1)-1
517     atom2 = groupList(jb)
518    
519     if (skipThisPair(atom1, atom2)) cycle prepair_inner
520    
521     call get_interatomic_vector(q_Row(:,atom1), &
522     q_Col(:,atom2), d_atm, ratmsq)
523    
524     call do_prepair(atom1, atom2, ratmsq, d_atm, &
525     rgrpsq, d_grp, do_pot, do_stress, &
526     u_l, A, f, t, pot_local)
527    
528     enddo prepair_inner
529     enddo
530     end if
531     enddo
532 gezelter 1138 enddo
533 gezelter 1150 point(nrow_group + 1) = nlist + 1
534    
535 gezelter 1138 else !! (of update_check)
536    
537     ! use the list to find the neighbors
538 gezelter 1150 do i = 1, nrow_group
539 gezelter 1138 JBEG = POINT(i)
540     JEND = POINT(i+1) - 1
541 gezelter 1150 ! check that group i has neighbors
542 gezelter 1138 if (jbeg .le. jend) then
543    
544     do jnab = jbeg, jend
545     j = list(jnab)
546    
547 gezelter 1150 do ia = groupStart(i), groupStart(i+1)-1
548     atom1 = groupList(ia)
549    
550     do jb = groupStart(j), groupStart(j+1)-1
551     atom2 = groupList(jb)
552    
553     call get_interatomic_vector(q_Row(:,atom1), &
554     q_Col(:,atom2), d_atm, ratmsq)
555    
556     call do_prepair(atom1, atom2, ratmsq, d_atm, &
557     rgrpsq, d_grp, do_pot, do_stress, &
558     u_l, A, f, t, pot_local)
559    
560     enddo
561     enddo
562 gezelter 1138 enddo
563 chuckv 648 endif
564 gezelter 1138 enddo
565     endif
566 chuckv 648
567     #else
568 gezelter 1150
569 gezelter 1138 if (update_nlist) then
570 chuckv 648
571 gezelter 1150 !! save current configuration, construct neighbor list,
572     !! and calculate forces
573 gezelter 1138
574 gezelter 1150 call saveNeighborList(nGroup, q_group)
575    
576 gezelter 1138 neighborListSize = size(list)
577 gezelter 1150 nlist = 0
578 gezelter 1138
579 gezelter 1150 do i = 1, nGroup-1
580 gezelter 1138 point(i) = nlist + 1
581 chuckv 648
582 gezelter 1150 do j = i+1, nGroup
583 gezelter 1138
584 gezelter 1150 call get_interatomic_vector(q_group(:,i), &
585     q_group(:,j), d_grp, rgrpsq)
586 gezelter 1138
587 gezelter 1150 if (rgrpsq < rlistsq) then
588 gezelter 1138 nlist = nlist + 1
589    
590     if (nlist > neighborListSize) then
591 gezelter 1150 call expandNeighborList(nGroup, listerror)
592 gezelter 1138 if (listerror /= 0) then
593     error = -1
594     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
595     return
596     end if
597     neighborListSize = size(list)
598     endif
599    
600     list(nlist) = j
601    
602 gezelter 1150 do ia = groupStart(i), groupStart(i+1)-1
603     atom1 = groupList(ia)
604    
605     prepair_inner: do jb = groupStart(j), groupStart(j+1)-1
606     atom2 = groupList(jb)
607    
608     if (skipThisPair(atom1, atom2)) cycle prepair_inner
609    
610     call get_interatomic_vector(q(:,atom1), &
611     q(:,atom2), d_atm, ratmsq)
612    
613     call do_prepair(atom1, atom2, ratmsq, d_atm, &
614     rgrpsq, d_grp, do_pot, do_stress, &
615     u_l, A, f, t, pot)
616    
617     enddo prepair_inner
618     enddo
619     end if
620     enddo
621 gezelter 1138 enddo
622 gezelter 1150 point(nGroup) = nlist + 1
623    
624     else !! (of update_check)
625 gezelter 1138
626     ! use the list to find the neighbors
627 gezelter 1150 do i = 1, nGroup-1
628 gezelter 1138 JBEG = POINT(i)
629     JEND = POINT(i+1) - 1
630 gezelter 1150 ! check that group i has neighbors
631 gezelter 1138 if (jbeg .le. jend) then
632 chuckv 648
633 gezelter 1138 do jnab = jbeg, jend
634     j = list(jnab)
635    
636 gezelter 1150 do ia = groupStart(i), groupStart(i+1)-1
637     atom1 = groupList(ia)
638    
639     do jb = groupStart(j), groupStart(j+1)-1
640     atom2 = groupList(jb)
641    
642     call get_interatomic_vector(q(:,atom1), &
643     q(:,atom2), d_atm, ratmsq)
644    
645     call do_prepair(atom1, atom2, ratmsq, d_atm, &
646     rgrpsq, d_grp, do_pot, do_stress, &
647     u_l, A, f, t, pot)
648    
649     enddo
650     enddo
651 gezelter 1138 enddo
652 chuckv 648 endif
653 gezelter 1138 enddo
654     endif
655 gezelter 1150
656 chuckv 648 #endif
657 gezelter 1150
658 gezelter 1138 !! Do rest of preforce calculations
659     !! do necessary preforce calculations
660     call do_preforce(nlocal,pot)
661     ! we have already updated the neighbor list set it to false...
662     update_nlist = .false.
663 mmeineke 377 else
664 chuckv 648 !! See if we need to update neighbor lists for non pre-pair
665 gezelter 1150 call checkNeighborList(nGroup, q_group, listSkin, update_nlist)
666 mmeineke 377 endif
667 gezelter 1138
668     !---------------------------------MAIN Pair LOOP->>>>>>>>>>>>>
669    
670 mmeineke 377 #ifdef IS_MPI
671    
672     if (update_nlist) then
673 gezelter 1150
674 mmeineke 377 !! save current configuration, construct neighbor list,
675     !! and calculate forces
676    
677 gezelter 1150 call saveNeighborList(nGroup, q_group)
678    
679 mmeineke 377 neighborListSize = size(list)
680     nlist = 0
681    
682 gezelter 1150 do i = 1, nrow_group
683 gezelter 1169
684 mmeineke 377 point(i) = nlist + 1
685 gezelter 1169
686     n_in_i = groupStart(i+1) - groupStart(i)
687 mmeineke 377
688 gezelter 1150 do j = 1, ncol_group
689 mmeineke 377
690 gezelter 1150 call get_interatomic_vector(q_group_Row(:,i), &
691     q_group_Col(:,j), d_grp, rgrpsq)
692 mmeineke 377
693 gezelter 1150 if (rgrpsq < rlistsq) then
694 mmeineke 377 nlist = nlist + 1
695    
696     if (nlist > neighborListSize) then
697 gezelter 1150 call expandNeighborList(nGroup, listerror)
698 mmeineke 377 if (listerror /= 0) then
699     error = -1
700     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
701     return
702     end if
703     neighborListSize = size(list)
704     endif
705    
706     list(nlist) = j
707 mmeineke 626
708 gezelter 1169 vij = 0.0d0
709    
710 gezelter 1150 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
711     in_switching_region)
712 gezelter 1169
713     n_in_j = groupStart(j+1) - groupStart(j)
714 gezelter 1150
715     do ia = groupStart(i), groupStart(i+1)-1
716     atom1 = groupList(ia)
717    
718     inner: do jb = groupStart(j), groupStart(j+1)-1
719 gezelter 1163 atom2 = groupList(jb)
720    
721 gezelter 1150 if (skipThisPair(atom1, atom2)) cycle inner
722    
723 gezelter 1169 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
724     d_atm(1:3) = d_grp(1:3)
725     ratmsq = rgrpsq
726     else
727     call get_interatomic_vector(q_Row(:,atom1), &
728     q_Col(:,atom2), d_atm, ratmsq)
729     endif
730 gezelter 1150
731     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
732     do_pot, do_stress, &
733     u_l, A, f, t, pot_local, vpair)
734    
735 gezelter 1169 vij = vij + vpair
736 gezelter 1150 enddo inner
737     enddo
738    
739     if (in_switching_region) then
740 gezelter 1169 swderiv = vij*dswdr/rgrp
741 gezelter 1150
742     do ia=groupStart(i), groupStart(i+1)-1
743     atom1=groupList(ia)
744     mf = mfact(atom1)
745 gezelter 1169 f_Row(1,atom1) = f_Row(1,atom1) + swderiv*d_grp(1)*mf
746     f_Row(2,atom1) = f_Row(2,atom1) + swderiv*d_grp(2)*mf
747     f_Row(3,atom1) = f_Row(3,atom1) + swderiv*d_grp(3)*mf
748 gezelter 1150 enddo
749    
750     do jb=groupStart(j), groupStart(j+1)-1
751     atom2=groupList(jb)
752     mf = mfact(atom2)
753 gezelter 1169 f_Col(1,atom2) = f_Col(1,atom2) - swderiv*d_grp(1)*mf
754     f_Col(2,atom2) = f_Col(2,atom2) - swderiv*d_grp(2)*mf
755     f_Col(3,atom2) = f_Col(3,atom2) - swderiv*d_grp(3)*mf
756 gezelter 1150 enddo
757 gezelter 1138 endif
758 gezelter 1150
759     end if
760     enddo
761 mmeineke 377 enddo
762 gezelter 1169
763 gezelter 1150 point(nrow_group + 1) = nlist + 1
764 gezelter 1138
765 mmeineke 377 else !! (of update_check)
766 gezelter 1138
767 mmeineke 377 ! use the list to find the neighbors
768 gezelter 1150 do i = 1, nrow_group
769 gezelter 1169
770     n_in_i = groupStart(i+1) - groupStart(i)
771    
772 mmeineke 377 JBEG = POINT(i)
773     JEND = POINT(i+1) - 1
774 gezelter 1150 ! check that group i has neighbors
775 mmeineke 377 if (jbeg .le. jend) then
776    
777     do jnab = jbeg, jend
778     j = list(jnab)
779 gezelter 1150
780     call get_interatomic_vector(q_group_Row(:,i), &
781     q_group_Col(:,j), d_grp, rgrpsq)
782 gezelter 1138
783 gezelter 1169 vij = 0.0d0
784 gezelter 1150 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
785     in_switching_region)
786 gezelter 1138
787 gezelter 1169 n_in_j = groupStart(j+1) - groupStart(j)
788    
789 gezelter 1150 do ia = groupStart(i), groupStart(i+1)-1
790     atom1 = groupList(ia)
791    
792     do jb = groupStart(j), groupStart(j+1)-1
793     atom2 = groupList(jb)
794    
795 gezelter 1169 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
796     d_atm(1:3) = d_grp(1:3)
797     ratmsq = rgrpsq
798     else
799     call get_interatomic_vector(q_Row(:,atom1), &
800     q_Col(:,atom2), d_atm, ratmsq)
801     endif
802 gezelter 1150
803     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
804     do_pot, do_stress, &
805     u_l, A, f, t, pot_local, vpair)
806    
807 gezelter 1169 vij = vij + vpair
808 gezelter 1150
809     enddo
810     enddo
811 gezelter 1138
812 gezelter 1150 if (in_switching_region) then
813 gezelter 1169 swderiv = vij*dswdr/rgrp
814 gezelter 1150
815     do ia=groupStart(i), groupStart(i+1)-1
816     atom1=groupList(ia)
817     mf = mfact(atom1)
818 gezelter 1169 f_Row(1,atom1) = f_Row(1,atom1) + swderiv*d_grp(1)*mf
819     f_Row(2,atom1) = f_Row(2,atom1) + swderiv*d_grp(2)*mf
820     f_Row(3,atom1) = f_Row(3,atom1) + swderiv*d_grp(3)*mf
821 gezelter 1150 enddo
822    
823     do jb=groupStart(j), groupStart(j+1)-1
824     atom2=groupList(jb)
825     mf = mfact(atom2)
826 gezelter 1169 f_Col(1,atom2) = f_Col(1,atom2) - swderiv*d_grp(1)*mf
827     f_Col(2,atom2) = f_Col(2,atom2) - swderiv*d_grp(2)*mf
828     f_Col(3,atom2) = f_Col(3,atom2) - swderiv*d_grp(3)*mf
829 gezelter 1150 enddo
830     endif
831    
832 mmeineke 377 enddo
833     endif
834     enddo
835     endif
836    
837     #else
838 gezelter 1150
839 mmeineke 377 if (update_nlist) then
840 gezelter 1138
841 gezelter 1150 !! save current configuration, construct neighbor list,
842     !! and calculate forces
843 mmeineke 377
844 gezelter 1150 call saveNeighborList(nGroup, q_group)
845    
846 mmeineke 377 neighborListSize = size(list)
847 gezelter 1150 nlist = 0
848 gezelter 1138
849 gezelter 1150 do i = 1, nGroup-1
850 gezelter 1169
851 mmeineke 377 point(i) = nlist + 1
852 gezelter 1169 n_in_i = groupStart(i+1) - groupStart(i)
853 mmeineke 377
854 gezelter 1150 do j = i+1, nGroup
855 mmeineke 377
856 gezelter 1150 call get_interatomic_vector(q_group(:,i), &
857     q_group(:,j), d_grp, rgrpsq)
858 gezelter 1138
859 gezelter 1150 if (rgrpsq < rlistsq) then
860 mmeineke 377 nlist = nlist + 1
861 gezelter 1138
862 mmeineke 377 if (nlist > neighborListSize) then
863 gezelter 1150 call expandNeighborList(nGroup, listerror)
864 mmeineke 377 if (listerror /= 0) then
865     error = -1
866     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
867     return
868     end if
869     neighborListSize = size(list)
870     endif
871    
872     list(nlist) = j
873 gezelter 1150
874 gezelter 1169 vij = 0.0d0
875    
876 gezelter 1150 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
877     in_switching_region)
878 mmeineke 377
879 gezelter 1169 n_in_j = groupStart(j+1) - groupStart(j)
880    
881 gezelter 1150 do ia = groupStart(i), groupStart(i+1)-1
882     atom1 = groupList(ia)
883    
884     inner: do jb = groupStart(j), groupStart(j+1)-1
885     atom2 = groupList(jb)
886    
887     if (skipThisPair(atom1, atom2)) cycle inner
888    
889 gezelter 1169 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
890     d_atm(1:3) = d_grp(1:3)
891     ratmsq = rgrpsq
892     else
893     call get_interatomic_vector(q(:,atom1), &
894     q(:,atom2), d_atm, ratmsq)
895     endif
896    
897 gezelter 1150 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
898     do_pot, do_stress, &
899     u_l, A, f, t, pot, vpair)
900    
901 gezelter 1169 vij = vij + vpair
902 gezelter 1150
903     enddo inner
904     enddo
905    
906     if (in_switching_region) then
907 gezelter 1169 swderiv = vij*dswdr/rgrp
908 gezelter 1150 do ia=groupStart(i), groupStart(i+1)-1
909     atom1=groupList(ia)
910     mf = mfact(atom1)
911 gezelter 1169 f(1,atom1) = f(1,atom1) + swderiv*d_grp(1)*mf
912     f(2,atom1) = f(2,atom1) + swderiv*d_grp(2)*mf
913     f(3,atom1) = f(3,atom1) + swderiv*d_grp(3)*mf
914 gezelter 1150 enddo
915    
916     do jb=groupStart(j), groupStart(j+1)-1
917     atom2=groupList(jb)
918     mf = mfact(atom2)
919 gezelter 1169 f(1,atom2) = f(1,atom2) - swderiv*d_grp(1)*mf
920     f(2,atom2) = f(2,atom2) - swderiv*d_grp(2)*mf
921     f(3,atom2) = f(3,atom2) - swderiv*d_grp(3)*mf
922 gezelter 1150 enddo
923 gezelter 1138 endif
924 gezelter 1150
925     end if
926     enddo
927 mmeineke 377 enddo
928 gezelter 1150 point(nGroup) = nlist + 1
929 mmeineke 377
930 gezelter 1150 else !! (of update_check)
931 mmeineke 377
932     ! use the list to find the neighbors
933 gezelter 1150 do i = 1, nGroup-1
934 gezelter 1169
935     n_in_i = groupStart(i+1) - groupStart(i)
936    
937 mmeineke 377 JBEG = POINT(i)
938     JEND = POINT(i+1) - 1
939 gezelter 1150 ! check that group i has neighbors
940 mmeineke 377 if (jbeg .le. jend) then
941    
942     do jnab = jbeg, jend
943     j = list(jnab)
944 gezelter 1169
945 gezelter 1150 call get_interatomic_vector(q_group(:,i), &
946     q_group(:,j), d_grp, rgrpsq)
947 gezelter 1138
948 gezelter 1169 vij = 0.0d0
949    
950 gezelter 1150 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
951     in_switching_region)
952 gezelter 1169
953     n_in_j = groupStart(j+1) - groupStart(j)
954    
955 gezelter 1150 do ia = groupStart(i), groupStart(i+1)-1
956     atom1 = groupList(ia)
957    
958     do jb = groupStart(j), groupStart(j+1)-1
959     atom2 = groupList(jb)
960    
961 gezelter 1169 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
962     d_atm(1:3) = d_grp(1:3)
963     ratmsq = rgrpsq
964     else
965     call get_interatomic_vector(q(:,atom1), &
966     q(:,atom2), d_atm, ratmsq)
967     endif
968 gezelter 1150
969     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
970     do_pot, do_stress, &
971     u_l, A, f, t, pot, vpair)
972    
973 gezelter 1169 vij = vij + vpair
974    
975 gezelter 1150 enddo
976     enddo
977    
978     if (in_switching_region) then
979 gezelter 1169 swderiv = vij*dswdr/rgrp
980 gezelter 1150
981     do ia=groupStart(i), groupStart(i+1)-1
982     atom1=groupList(ia)
983     mf = mfact(atom1)
984 gezelter 1169 f(1,atom1) = f(1,atom1) + swderiv*d_grp(1)*mf
985     f(2,atom1) = f(2,atom1) + swderiv*d_grp(2)*mf
986     f(3,atom1) = f(3,atom1) + swderiv*d_grp(3)*mf
987 gezelter 1150 enddo
988    
989     do jb=groupStart(j), groupStart(j+1)-1
990     atom2=groupList(jb)
991     mf = mfact(atom2)
992 gezelter 1169 f(1,atom2) = f(1,atom2) - swderiv*d_grp(1)*mf
993     f(2,atom2) = f(2,atom2) - swderiv*d_grp(2)*mf
994     f(3,atom2) = f(3,atom2) - swderiv*d_grp(3)*mf
995 gezelter 1150 enddo
996 gezelter 1138 endif
997 mmeineke 377 enddo
998     endif
999     enddo
1000     endif
1001 gezelter 1169
1002 mmeineke 377 #endif
1003    
1004     ! phew, done with main loop.
1005 gezelter 1138
1006     !! Do timing
1007 chuckv 694 #ifdef PROFILE
1008     call cpu_time(forceTimeFinal)
1009     forceTime = forceTime + forceTimeFinal - forceTimeInitial
1010 gezelter 1150 #endif
1011 gezelter 1138
1012 mmeineke 377 #ifdef IS_MPI
1013     !!distribute forces
1014 gezelter 1138
1015 chuckv 438 f_temp = 0.0_dp
1016     call scatter(f_Row,f_temp,plan_row3d)
1017     do i = 1,nlocal
1018     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1019     end do
1020 gezelter 1138
1021 chuckv 438 f_temp = 0.0_dp
1022 mmeineke 377 call scatter(f_Col,f_temp,plan_col3d)
1023     do i = 1,nlocal
1024     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
1025     end do
1026    
1027 chuckv 900 if (FF_UsesDirectionalAtoms() .and. SIM_uses_directional_atoms) then
1028 chuckv 438 t_temp = 0.0_dp
1029     call scatter(t_Row,t_temp,plan_row3d)
1030     do i = 1,nlocal
1031     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1032     end do
1033     t_temp = 0.0_dp
1034 mmeineke 377 call scatter(t_Col,t_temp,plan_col3d)
1035    
1036     do i = 1,nlocal
1037     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
1038     end do
1039     endif
1040    
1041     if (do_pot) then
1042     ! scatter/gather pot_row into the members of my column
1043     call scatter(pot_Row, pot_Temp, plan_row)
1044 gezelter 1138
1045 mmeineke 377 ! scatter/gather pot_local into all other procs
1046     ! add resultant to get total pot
1047     do i = 1, nlocal
1048     pot_local = pot_local + pot_Temp(i)
1049     enddo
1050 chuckv 439
1051     pot_Temp = 0.0_DP
1052 gezelter 1138
1053 mmeineke 377 call scatter(pot_Col, pot_Temp, plan_col)
1054     do i = 1, nlocal
1055     pot_local = pot_local + pot_Temp(i)
1056     enddo
1057 gezelter 1150
1058 gezelter 1138 endif
1059 mmeineke 377 #endif
1060 gezelter 1138
1061 chuckv 900 if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then
1062 mmeineke 377
1063 chuckv 900 if (FF_uses_RF .and. SIM_uses_RF) then
1064 gezelter 1138
1065 mmeineke 377 #ifdef IS_MPI
1066     call scatter(rf_Row,rf,plan_row3d)
1067     call scatter(rf_Col,rf_Temp,plan_col3d)
1068     do i = 1,nlocal
1069     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
1070     end do
1071     #endif
1072    
1073 chuckv 898 do i = 1, nLocal
1074 gezelter 1138
1075 mmeineke 377 rfpot = 0.0_DP
1076     #ifdef IS_MPI
1077     me_i = atid_row(i)
1078     #else
1079     me_i = atid(i)
1080     #endif
1081 gezelter 1138
1082 chuckv 900 if (PropertyMap(me_i)%is_DP) then
1083 gezelter 1138
1084 chuckv 900 mu_i = PropertyMap(me_i)%dipole_moment
1085 gezelter 1138
1086 mmeineke 377 !! The reaction field needs to include a self contribution
1087     !! to the field:
1088 chuckv 900 call accumulate_self_rf(i, mu_i, u_l)
1089 mmeineke 377 !! Get the reaction field contribution to the
1090     !! potential and torques:
1091     call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
1092     #ifdef IS_MPI
1093     pot_local = pot_local + rfpot
1094     #else
1095     pot = pot + rfpot
1096    
1097     #endif
1098     endif
1099     enddo
1100     endif
1101     endif
1102 gezelter 1138
1103    
1104 mmeineke 377 #ifdef IS_MPI
1105 gezelter 1138
1106 mmeineke 377 if (do_pot) then
1107 chuckv 441 pot = pot + pot_local
1108 mmeineke 377 !! we assume the c code will do the allreduce to get the total potential
1109     !! we could do it right here if we needed to...
1110     endif
1111 gezelter 1138
1112 mmeineke 377 if (do_stress) then
1113 gezelter 1138 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
1114 mmeineke 377 mpi_comm_world,mpi_err)
1115 chuckv 470 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
1116 mmeineke 377 mpi_comm_world,mpi_err)
1117     endif
1118 gezelter 1138
1119 mmeineke 377 #else
1120 gezelter 1138
1121 mmeineke 377 if (do_stress) then
1122     tau = tau_Temp
1123     virial = virial_Temp
1124     endif
1125 mmeineke 887
1126 mmeineke 377 #endif
1127 mmeineke 887
1128    
1129 mmeineke 377 end subroutine do_force_loop
1130 gezelter 1150
1131 gezelter 1138
1132 gezelter 1150 subroutine do_pair(i, j, rijsq, d, sw, do_pot, do_stress, &
1133     u_l, A, f, t, pot, vpair)
1134 mmeineke 377
1135 gezelter 1150 real( kind = dp ) :: pot, vpair, sw
1136 gezelter 1138 real( kind = dp ), dimension(nLocal) :: mfact
1137 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
1138 gezelter 1138 real( kind = dp ), dimension(9,nLocal) :: A
1139     real( kind = dp ), dimension(3,nLocal) :: f
1140     real( kind = dp ), dimension(3,nLocal) :: t
1141 mmeineke 377
1142     logical, intent(inout) :: do_pot, do_stress
1143     integer, intent(in) :: i, j
1144 gezelter 1150 real ( kind = dp ), intent(inout) :: rijsq
1145     real ( kind = dp ) :: r
1146     real ( kind = dp ), intent(inout) :: d(3)
1147 mmeineke 377 integer :: me_i, me_j
1148 gezelter 946
1149 mmeineke 377 r = sqrt(rijsq)
1150 gezelter 1163 vpair = 0.0d0
1151 mmeineke 377
1152     #ifdef IS_MPI
1153 gezelter 490 if (tagRow(i) .eq. tagColumn(j)) then
1154     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
1155     endif
1156 mmeineke 377 me_i = atid_row(i)
1157     me_j = atid_col(j)
1158     #else
1159     me_i = atid(i)
1160     me_j = atid(j)
1161     #endif
1162 chuckv 894
1163 chuckv 900 if (FF_uses_LJ .and. SIM_uses_LJ) then
1164 gezelter 946
1165     if ( PropertyMap(me_i)%is_LJ .and. PropertyMap(me_j)%is_LJ ) then
1166 gezelter 1163 !write(*,*) 'calling lj with'
1167     !write(*,*) i, j, r, rijsq
1168     !write(*,'(3es12.3)') d(1), d(2), d(3)
1169     !write(*,'(3es12.3)') sw, vpair, pot
1170     !write(*,*)
1171    
1172 gezelter 1150 call do_lj_pair(i, j, d, r, rijsq, sw, vpair, pot, f, do_pot, &
1173     do_stress)
1174 gezelter 946 endif
1175    
1176 mmeineke 377 endif
1177 gezelter 946
1178     if (FF_uses_charges .and. SIM_uses_charges) then
1179    
1180     if (PropertyMap(me_i)%is_Charge .and. PropertyMap(me_j)%is_Charge) then
1181 gezelter 1150 call do_charge_pair(i, j, d, r, rijsq, sw, vpair, pot, f, do_pot, &
1182     do_stress)
1183 gezelter 946 endif
1184    
1185     endif
1186    
1187 chuckv 900 if (FF_uses_dipoles .and. SIM_uses_dipoles) then
1188 mmeineke 377
1189 chuckv 900 if ( PropertyMap(me_i)%is_DP .and. PropertyMap(me_j)%is_DP) then
1190 gezelter 1150 call do_dipole_pair(i, j, d, r, rijsq, sw, vpair, pot, u_l, f, t, &
1191 mmeineke 377 do_pot, do_stress)
1192 chuckv 900 if (FF_uses_RF .and. SIM_uses_RF) then
1193 gezelter 1160 call accumulate_rf(i, j, r, u_l, sw)
1194     call rf_correct_forces(i, j, d, r, u_l, sw, f, do_stress)
1195 gezelter 946 endif
1196 mmeineke 377 endif
1197 gezelter 946
1198 mmeineke 377 endif
1199    
1200 chuckv 900 if (FF_uses_Sticky .and. SIM_uses_sticky) then
1201 mmeineke 377
1202 chuckv 900 if ( PropertyMap(me_i)%is_Sticky .and. PropertyMap(me_j)%is_Sticky) then
1203 gezelter 1150 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, pot, A, f, t, &
1204 mmeineke 377 do_pot, do_stress)
1205     endif
1206 gezelter 946
1207 mmeineke 377 endif
1208    
1209    
1210 chuckv 900 if (FF_uses_GB .and. SIM_uses_GB) then
1211 mmeineke 377
1212 chuckv 900 if ( PropertyMap(me_i)%is_GB .and. PropertyMap(me_j)%is_GB) then
1213 gezelter 1150 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, pot, u_l, f, t, &
1214 mmeineke 377 do_pot, do_stress)
1215     endif
1216 chuckv 894
1217 mmeineke 377 endif
1218 gezelter 946
1219     if (FF_uses_EAM .and. SIM_uses_EAM) then
1220    
1221     if ( PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) then
1222 gezelter 1150 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, pot, f, &
1223     do_pot, do_stress)
1224 gezelter 946 endif
1225    
1226     endif
1227 gezelter 1150
1228 mmeineke 377 end subroutine do_pair
1229    
1230 gezelter 1138 subroutine do_prepair(i, j, rijsq, d, rcijsq, dc, &
1231     do_pot, do_stress, u_l, A, f, t, pot)
1232 chuckv 631 real( kind = dp ) :: pot
1233 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
1234     real (kind=dp), dimension(9,nLocal) :: A
1235     real (kind=dp), dimension(3,nLocal) :: f
1236     real (kind=dp), dimension(3,nLocal) :: t
1237 chuckv 631
1238     logical, intent(inout) :: do_pot, do_stress
1239     integer, intent(in) :: i, j
1240 gezelter 1138 real ( kind = dp ), intent(inout) :: rijsq, rcijsq
1241     real ( kind = dp ) :: r, rc
1242     real ( kind = dp ), intent(inout) :: d(3), dc(3)
1243 chuckv 631
1244     logical :: is_EAM_i, is_EAM_j
1245    
1246     integer :: me_i, me_j
1247    
1248 gezelter 1138
1249     r = sqrt(rijsq)
1250     if (SIM_uses_molecular_cutoffs) then
1251     rc = sqrt(rcijsq)
1252     else
1253     rc = r
1254     endif
1255 chuckv 631
1256 chuckv 669
1257 chuckv 631 #ifdef IS_MPI
1258     if (tagRow(i) .eq. tagColumn(j)) then
1259     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
1260     endif
1261    
1262     me_i = atid_row(i)
1263     me_j = atid_col(j)
1264    
1265     #else
1266    
1267     me_i = atid(i)
1268     me_j = atid(j)
1269    
1270     #endif
1271 chuckv 673
1272 chuckv 900 if (FF_uses_EAM .and. SIM_uses_EAM) then
1273    
1274     if (PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) &
1275 chuckv 648 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
1276 chuckv 900
1277 chuckv 631 endif
1278 chuckv 900
1279 chuckv 673 end subroutine do_prepair
1280 chuckv 648
1281    
1282    
1283 chuckv 673
1284 gezelter 1138 subroutine do_preforce(nlocal,pot)
1285     integer :: nlocal
1286     real( kind = dp ) :: pot
1287    
1288     if (FF_uses_EAM .and. SIM_uses_EAM) then
1289     call calc_EAM_preforce_Frho(nlocal,pot)
1290     endif
1291    
1292    
1293     end subroutine do_preforce
1294    
1295    
1296     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1297    
1298     real (kind = dp), dimension(3) :: q_i
1299     real (kind = dp), dimension(3) :: q_j
1300     real ( kind = dp ), intent(out) :: r_sq
1301     real( kind = dp ) :: d(3), scaled(3)
1302     integer i
1303    
1304     d(1:3) = q_j(1:3) - q_i(1:3)
1305    
1306     ! Wrap back into periodic box if necessary
1307     if ( SIM_uses_PBC ) then
1308 mmeineke 377
1309 gezelter 1138 if( .not.boxIsOrthorhombic ) then
1310     ! calc the scaled coordinates.
1311    
1312     scaled = matmul(HmatInv, d)
1313    
1314     ! wrap the scaled coordinates
1315    
1316     scaled = scaled - anint(scaled)
1317    
1318    
1319     ! calc the wrapped real coordinates from the wrapped scaled
1320     ! coordinates
1321    
1322     d = matmul(Hmat,scaled)
1323    
1324     else
1325     ! calc the scaled coordinates.
1326    
1327     do i = 1, 3
1328     scaled(i) = d(i) * HmatInv(i,i)
1329    
1330     ! wrap the scaled coordinates
1331    
1332     scaled(i) = scaled(i) - anint(scaled(i))
1333    
1334     ! calc the wrapped real coordinates from the wrapped scaled
1335     ! coordinates
1336    
1337     d(i) = scaled(i)*Hmat(i,i)
1338     enddo
1339     endif
1340    
1341     endif
1342    
1343     r_sq = dot_product(d,d)
1344    
1345     end subroutine get_interatomic_vector
1346 mmeineke 377
1347 gezelter 1138 subroutine zero_work_arrays()
1348    
1349     #ifdef IS_MPI
1350    
1351     q_Row = 0.0_dp
1352     q_Col = 0.0_dp
1353 mmeineke 377
1354 gezelter 1150 q_group_Row = 0.0_dp
1355     q_group_Col = 0.0_dp
1356 gezelter 1138
1357     u_l_Row = 0.0_dp
1358     u_l_Col = 0.0_dp
1359    
1360     A_Row = 0.0_dp
1361     A_Col = 0.0_dp
1362    
1363     f_Row = 0.0_dp
1364     f_Col = 0.0_dp
1365     f_Temp = 0.0_dp
1366    
1367     t_Row = 0.0_dp
1368     t_Col = 0.0_dp
1369     t_Temp = 0.0_dp
1370    
1371     pot_Row = 0.0_dp
1372     pot_Col = 0.0_dp
1373     pot_Temp = 0.0_dp
1374    
1375     rf_Row = 0.0_dp
1376     rf_Col = 0.0_dp
1377     rf_Temp = 0.0_dp
1378    
1379 mmeineke 377 #endif
1380 chuckv 673
1381 gezelter 1138 if (FF_uses_EAM .and. SIM_uses_EAM) then
1382     call clean_EAM()
1383     endif
1384    
1385     rf = 0.0_dp
1386     tau_Temp = 0.0_dp
1387     virial_Temp = 0.0_dp
1388     end subroutine zero_work_arrays
1389    
1390     function skipThisPair(atom1, atom2) result(skip_it)
1391     integer, intent(in) :: atom1
1392     integer, intent(in), optional :: atom2
1393     logical :: skip_it
1394     integer :: unique_id_1, unique_id_2
1395     integer :: me_i,me_j
1396     integer :: i
1397    
1398     skip_it = .false.
1399    
1400     !! there are a number of reasons to skip a pair or a particle
1401     !! mostly we do this to exclude atoms who are involved in short
1402     !! range interactions (bonds, bends, torsions), but we also need
1403     !! to exclude some overcounted interactions that result from
1404     !! the parallel decomposition
1405    
1406 mmeineke 377 #ifdef IS_MPI
1407 gezelter 1138 !! in MPI, we have to look up the unique IDs for each atom
1408     unique_id_1 = tagRow(atom1)
1409 mmeineke 377 #else
1410 gezelter 1138 !! in the normal loop, the atom numbers are unique
1411     unique_id_1 = atom1
1412 mmeineke 377 #endif
1413 gezelter 1138
1414     !! We were called with only one atom, so just check the global exclude
1415     !! list for this atom
1416     if (.not. present(atom2)) then
1417     do i = 1, nExcludes_global
1418     if (excludesGlobal(i) == unique_id_1) then
1419     skip_it = .true.
1420     return
1421     end if
1422     end do
1423     return
1424     end if
1425    
1426 mmeineke 377 #ifdef IS_MPI
1427 gezelter 1138 unique_id_2 = tagColumn(atom2)
1428 mmeineke 377 #else
1429 gezelter 1138 unique_id_2 = atom2
1430 mmeineke 377 #endif
1431 gezelter 1138
1432 mmeineke 377 #ifdef IS_MPI
1433 gezelter 1138 !! this situation should only arise in MPI simulations
1434     if (unique_id_1 == unique_id_2) then
1435     skip_it = .true.
1436     return
1437     end if
1438    
1439     !! this prevents us from doing the pair on multiple processors
1440     if (unique_id_1 < unique_id_2) then
1441     if (mod(unique_id_1 + unique_id_2,2) == 0) then
1442     skip_it = .true.
1443     return
1444     endif
1445     else
1446     if (mod(unique_id_1 + unique_id_2,2) == 1) then
1447     skip_it = .true.
1448     return
1449     endif
1450     endif
1451 mmeineke 377 #endif
1452 gezelter 1138
1453     !! the rest of these situations can happen in all simulations:
1454     do i = 1, nExcludes_global
1455     if ((excludesGlobal(i) == unique_id_1) .or. &
1456     (excludesGlobal(i) == unique_id_2)) then
1457     skip_it = .true.
1458     return
1459     endif
1460     enddo
1461    
1462     do i = 1, nExcludes_local
1463     if (excludesLocal(1,i) == unique_id_1) then
1464     if (excludesLocal(2,i) == unique_id_2) then
1465     skip_it = .true.
1466     return
1467     endif
1468     else
1469     if (excludesLocal(1,i) == unique_id_2) then
1470     if (excludesLocal(2,i) == unique_id_1) then
1471     skip_it = .true.
1472     return
1473     endif
1474     endif
1475     endif
1476     end do
1477    
1478     return
1479     end function skipThisPair
1480 chuckv 441
1481 gezelter 1138 function FF_UsesDirectionalAtoms() result(doesit)
1482     logical :: doesit
1483     doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
1484     FF_uses_GB .or. FF_uses_RF
1485     end function FF_UsesDirectionalAtoms
1486    
1487     function FF_RequiresPrepairCalc() result(doesit)
1488     logical :: doesit
1489     doesit = FF_uses_EAM
1490     end function FF_RequiresPrepairCalc
1491    
1492     function FF_RequiresPostpairCalc() result(doesit)
1493     logical :: doesit
1494     doesit = FF_uses_RF
1495     end function FF_RequiresPostpairCalc
1496    
1497 chuckv 883 #ifdef PROFILE
1498 gezelter 1138 function getforcetime() result(totalforcetime)
1499     real(kind=dp) :: totalforcetime
1500     totalforcetime = forcetime
1501     end function getforcetime
1502 chuckv 883 #endif
1503 gezelter 1138
1504     !! This cleans componets of force arrays belonging only to fortran
1505    
1506 mmeineke 377 end module do_Forces