ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 1217
Committed: Tue Jun 1 21:45:22 2004 UTC (20 years, 1 month ago) by gezelter
File size: 33777 byte(s)
Log Message:
Bug fix (fixes of skipList and neighbor list under MPI)

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 1217 !! @version $Id: do_Forces.F90,v 1.68 2004-06-01 21:45:22 gezelter Exp $, $Date: 2004-06-01 21:45:22 $, $Name: not supported by cvs2svn $, $Revision: 1.68 $
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 gezelter 1197 INTEGER, PARAMETER:: PREPAIR_LOOP = 1
37     INTEGER, PARAMETER:: PAIR_LOOP = 2
38    
39 chuckv 900 logical, save :: haveRlist = .false.
40     logical, save :: haveNeighborList = .false.
41 mmeineke 626 logical, save :: havePolicies = .false.
42 chuckv 900 logical, save :: haveSIMvariables = .false.
43     logical, save :: havePropertyMap = .false.
44     logical, save :: haveSaneForceField = .false.
45 mmeineke 377 logical, save :: FF_uses_LJ
46     logical, save :: FF_uses_sticky
47 gezelter 941 logical, save :: FF_uses_charges
48 mmeineke 377 logical, save :: FF_uses_dipoles
49     logical, save :: FF_uses_RF
50     logical, save :: FF_uses_GB
51     logical, save :: FF_uses_EAM
52 chuckv 900 logical, save :: SIM_uses_LJ
53     logical, save :: SIM_uses_sticky
54 gezelter 941 logical, save :: SIM_uses_charges
55 chuckv 900 logical, save :: SIM_uses_dipoles
56     logical, save :: SIM_uses_RF
57     logical, save :: SIM_uses_GB
58     logical, save :: SIM_uses_EAM
59     logical, save :: SIM_requires_postpair_calc
60     logical, save :: SIM_requires_prepair_calc
61     logical, save :: SIM_uses_directional_atoms
62     logical, save :: SIM_uses_PBC
63 gezelter 1138 logical, save :: SIM_uses_molecular_cutoffs
64 mmeineke 377
65 mmeineke 626 real(kind=dp), save :: rlist, rlistsq
66    
67 mmeineke 377 public :: init_FF
68     public :: do_force_loop
69 mmeineke 626 public :: setRlistDF
70 mmeineke 377
71 chuckv 694 #ifdef PROFILE
72 chuckv 883 public :: getforcetime
73     real, save :: forceTime = 0
74     real :: forceTimeInitial, forceTimeFinal
75 mmeineke 891 integer :: nLoops
76 chuckv 694 #endif
77    
78 chuckv 900 type :: Properties
79     logical :: is_lj = .false.
80     logical :: is_sticky = .false.
81     logical :: is_dp = .false.
82     logical :: is_gb = .false.
83     logical :: is_eam = .false.
84 gezelter 941 logical :: is_charge = .false.
85     real(kind=DP) :: charge = 0.0_DP
86 chuckv 900 real(kind=DP) :: dipole_moment = 0.0_DP
87     end type Properties
88 chuckv 895
89 chuckv 900 type(Properties), dimension(:),allocatable :: PropertyMap
90    
91 mmeineke 377 contains
92    
93 mmeineke 626 subroutine setRlistDF( this_rlist )
94    
95     real(kind=dp) :: this_rlist
96    
97     rlist = this_rlist
98     rlistsq = rlist * rlist
99    
100     haveRlist = .true.
101    
102     end subroutine setRlistDF
103    
104 chuckv 900 subroutine createPropertyMap(status)
105     integer :: nAtypes
106     integer :: status
107     integer :: i
108     logical :: thisProperty
109     real (kind=DP) :: thisDPproperty
110    
111     status = 0
112    
113     nAtypes = getSize(atypes)
114    
115     if (nAtypes == 0) then
116     status = -1
117     return
118     end if
119    
120     if (.not. allocated(PropertyMap)) then
121     allocate(PropertyMap(nAtypes))
122     endif
123    
124     do i = 1, nAtypes
125     call getElementProperty(atypes, i, "is_LJ", thisProperty)
126     PropertyMap(i)%is_LJ = thisProperty
127 gezelter 941
128     call getElementProperty(atypes, i, "is_Charge", thisProperty)
129     PropertyMap(i)%is_Charge = thisProperty
130    
131     if (thisProperty) then
132     call getElementProperty(atypes, i, "charge", thisDPproperty)
133     PropertyMap(i)%charge = thisDPproperty
134     endif
135    
136 chuckv 900 call getElementProperty(atypes, i, "is_DP", thisProperty)
137     PropertyMap(i)%is_DP = thisProperty
138    
139     if (thisProperty) then
140     call getElementProperty(atypes, i, "dipole_moment", thisDPproperty)
141     PropertyMap(i)%dipole_moment = thisDPproperty
142     endif
143    
144     call getElementProperty(atypes, i, "is_Sticky", thisProperty)
145     PropertyMap(i)%is_Sticky = thisProperty
146     call getElementProperty(atypes, i, "is_GB", thisProperty)
147     PropertyMap(i)%is_GB = thisProperty
148     call getElementProperty(atypes, i, "is_EAM", thisProperty)
149     PropertyMap(i)%is_EAM = thisProperty
150     end do
151    
152     havePropertyMap = .true.
153    
154     end subroutine createPropertyMap
155    
156     subroutine setSimVariables()
157     SIM_uses_LJ = SimUsesLJ()
158     SIM_uses_sticky = SimUsesSticky()
159 gezelter 941 SIM_uses_charges = SimUsesCharges()
160 chuckv 900 SIM_uses_dipoles = SimUsesDipoles()
161     SIM_uses_RF = SimUsesRF()
162     SIM_uses_GB = SimUsesGB()
163     SIM_uses_EAM = SimUsesEAM()
164     SIM_requires_postpair_calc = SimRequiresPostpairCalc()
165     SIM_requires_prepair_calc = SimRequiresPrepairCalc()
166     SIM_uses_directional_atoms = SimUsesDirectionalAtoms()
167     SIM_uses_PBC = SimUsesPBC()
168 tim 1144 !SIM_uses_molecular_cutoffs = SimUsesMolecularCutoffs()
169 chuckv 900
170     haveSIMvariables = .true.
171    
172     return
173     end subroutine setSimVariables
174    
175     subroutine doReadyCheck(error)
176     integer, intent(out) :: error
177    
178     integer :: myStatus
179    
180     error = 0
181    
182     if (.not. havePropertyMap) then
183    
184     myStatus = 0
185    
186     call createPropertyMap(myStatus)
187    
188     if (myStatus .ne. 0) then
189     write(default_error, *) 'createPropertyMap failed in do_Forces!'
190     error = -1
191     return
192     endif
193     endif
194    
195     if (.not. haveSIMvariables) then
196     call setSimVariables()
197     endif
198    
199     if (.not. haveRlist) then
200     write(default_error, *) 'rList has not been set in do_Forces!'
201     error = -1
202     return
203     endif
204    
205     if (SIM_uses_LJ .and. FF_uses_LJ) then
206     if (.not. havePolicies) then
207     write(default_error, *) 'LJ mixing Policies have not been set in do_Forces!'
208     error = -1
209     return
210     endif
211     endif
212    
213     if (.not. haveNeighborList) then
214     write(default_error, *) 'neighbor list has not been initialized in do_Forces!'
215     error = -1
216     return
217     end if
218    
219     if (.not. haveSaneForceField) then
220     write(default_error, *) 'Force Field is not sane in do_Forces!'
221     error = -1
222     return
223     end if
224    
225     #ifdef IS_MPI
226     if (.not. isMPISimSet()) then
227     write(default_error,*) "ERROR: mpiSimulation has not been initialized!"
228     error = -1
229     return
230     endif
231     #endif
232     return
233     end subroutine doReadyCheck
234    
235    
236 mmeineke 377 subroutine init_FF(LJMIXPOLICY, use_RF_c, thisStat)
237    
238     integer, intent(in) :: LJMIXPOLICY
239     logical, intent(in) :: use_RF_c
240    
241     integer, intent(out) :: thisStat
242     integer :: my_status, nMatches
243     integer, pointer :: MatchList(:) => null()
244     real(kind=dp) :: rcut, rrf, rt, dielect
245    
246     !! assume things are copacetic, unless they aren't
247     thisStat = 0
248    
249     !! Fortran's version of a cast:
250     FF_uses_RF = use_RF_c
251    
252     !! init_FF is called *after* all of the atom types have been
253     !! defined in atype_module using the new_atype subroutine.
254     !!
255     !! this will scan through the known atypes and figure out what
256     !! interactions are used by the force field.
257    
258     FF_uses_LJ = .false.
259     FF_uses_sticky = .false.
260 gezelter 941 FF_uses_charges = .false.
261 mmeineke 377 FF_uses_dipoles = .false.
262     FF_uses_GB = .false.
263     FF_uses_EAM = .false.
264    
265     call getMatchingElementList(atypes, "is_LJ", .true., nMatches, MatchList)
266     if (nMatches .gt. 0) FF_uses_LJ = .true.
267 gezelter 941
268     call getMatchingElementList(atypes, "is_Charge", .true., nMatches, MatchList)
269 tim 1113 if (nMatches .gt. 0) FF_uses_charges = .true.
270    
271 mmeineke 377 call getMatchingElementList(atypes, "is_DP", .true., nMatches, MatchList)
272     if (nMatches .gt. 0) FF_uses_dipoles = .true.
273    
274     call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, &
275     MatchList)
276     if (nMatches .gt. 0) FF_uses_Sticky = .true.
277    
278     call getMatchingElementList(atypes, "is_GB", .true., nMatches, MatchList)
279     if (nMatches .gt. 0) FF_uses_GB = .true.
280    
281     call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList)
282     if (nMatches .gt. 0) FF_uses_EAM = .true.
283    
284 chuckv 900 !! Assume sanity (for the sake of argument)
285     haveSaneForceField = .true.
286    
287 mmeineke 377 !! check to make sure the FF_uses_RF setting makes sense
288    
289     if (FF_uses_dipoles) then
290     if (FF_uses_RF) then
291     dielect = getDielect()
292 mmeineke 626 call initialize_rf(dielect)
293 mmeineke 377 endif
294     else
295     if (FF_uses_RF) then
296     write(default_error,*) 'Using Reaction Field with no dipoles? Huh?'
297     thisStat = -1
298 chuckv 900 haveSaneForceField = .false.
299 mmeineke 377 return
300     endif
301 mmeineke 626 endif
302 mmeineke 377
303     if (FF_uses_LJ) then
304    
305     select case (LJMIXPOLICY)
306 gezelter 834 case (LB_MIXING_RULE)
307     call init_lj_FF(LB_MIXING_RULE, my_status)
308     case (EXPLICIT_MIXING_RULE)
309     call init_lj_FF(EXPLICIT_MIXING_RULE, my_status)
310 mmeineke 377 case default
311     write(default_error,*) 'unknown LJ Mixing Policy!'
312     thisStat = -1
313 chuckv 900 haveSaneForceField = .false.
314 mmeineke 377 return
315     end select
316     if (my_status /= 0) then
317     thisStat = -1
318 chuckv 900 haveSaneForceField = .false.
319 mmeineke 377 return
320     end if
321 chuckv 900 havePolicies = .true.
322 mmeineke 377 endif
323    
324     if (FF_uses_sticky) then
325     call check_sticky_FF(my_status)
326     if (my_status /= 0) then
327     thisStat = -1
328 chuckv 900 haveSaneForceField = .false.
329 mmeineke 377 return
330     end if
331     endif
332 chuckv 657
333    
334     if (FF_uses_EAM) then
335 chuckv 801 call init_EAM_FF(my_status)
336 chuckv 657 if (my_status /= 0) then
337 chuckv 900 write(default_error, *) "init_EAM_FF returned a bad status"
338 chuckv 657 thisStat = -1
339 chuckv 900 haveSaneForceField = .false.
340 chuckv 657 return
341     end if
342     endif
343    
344 mmeineke 377 if (FF_uses_GB) then
345     call check_gb_pair_FF(my_status)
346     if (my_status .ne. 0) then
347     thisStat = -1
348 chuckv 900 haveSaneForceField = .false.
349 mmeineke 377 return
350     endif
351     endif
352    
353     if (FF_uses_GB .and. FF_uses_LJ) then
354     endif
355 chuckv 900 if (.not. haveNeighborList) then
356 chuckv 480 !! Create neighbor lists
357 chuckv 898 call expandNeighborList(nLocal, my_status)
358 chuckv 480 if (my_Status /= 0) then
359     write(default_error,*) "SimSetup: ExpandNeighborList returned error."
360     thisStat = -1
361     return
362     endif
363 chuckv 900 haveNeighborList = .true.
364 chuckv 480 endif
365 gezelter 1150
366 chuckv 657
367 gezelter 1150
368 mmeineke 377 end subroutine init_FF
369    
370    
371     !! Does force loop over i,j pairs. Calls do_pair to calculates forces.
372     !------------------------------------------------------------->
373 gezelter 1150 subroutine do_force_loop(q, q_group, A, u_l, f, t, tau, pot, &
374 gezelter 1138 do_pot_c, do_stress_c, error)
375 mmeineke 377 !! Position array provided by C, dimensioned by getNlocal
376 gezelter 1166 real ( kind = dp ), dimension(3, nLocal) :: q
377 gezelter 1138 !! molecular center-of-mass position array
378 tim 1198 real ( kind = dp ), dimension(3, nGroups) :: q_group
379 mmeineke 377 !! Rotation Matrix for each long range particle in simulation.
380 gezelter 1166 real( kind = dp), dimension(9, nLocal) :: A
381 mmeineke 377 !! Unit vectors for dipoles (lab frame)
382 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
383 mmeineke 377 !! Force array provided by C, dimensioned by getNlocal
384 chuckv 898 real ( kind = dp ), dimension(3,nLocal) :: f
385 mmeineke 377 !! Torsion array provided by C, dimensioned by getNlocal
386 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: t
387 chuckv 895
388 mmeineke 377 !! Stress Tensor
389     real( kind = dp), dimension(9) :: tau
390     real ( kind = dp ) :: pot
391     logical ( kind = 2) :: do_pot_c, do_stress_c
392     logical :: do_pot
393     logical :: do_stress
394 gezelter 1150 logical :: in_switching_region
395 chuckv 439 #ifdef IS_MPI
396 chuckv 441 real( kind = DP ) :: pot_local
397 tim 1198 integer :: nAtomsInRow
398     integer :: nAtomsInCol
399 chuckv 694 integer :: nprocs
400 tim 1198 integer :: nGroupsInRow
401     integer :: nGroupsInCol
402 mmeineke 377 #endif
403     integer :: natoms
404     logical :: update_nlist
405 gezelter 1197 integer :: i, j, jstart, jend, jnab
406     integer :: istart, iend
407 gezelter 1150 integer :: ia, jb, atom1, atom2
408 mmeineke 377 integer :: nlist
409 gezelter 1169 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, vpair, vij
410 gezelter 1150 real( kind = DP ) :: sw, dswdr, swderiv, mf
411 gezelter 1192 real(kind=dp),dimension(3) :: d_atm, d_grp, fpair, fij
412 mmeineke 377 real(kind=dp) :: rfpot, mu_i, virial
413 gezelter 1169 integer :: me_i, me_j, n_in_i, n_in_j
414 mmeineke 377 logical :: is_dp_i
415     integer :: neighborListSize
416     integer :: listerror, error
417     integer :: localError
418 chuckv 897 integer :: propPack_i, propPack_j
419 gezelter 1197 integer :: loopStart, loopEnd, loop
420 mmeineke 626
421 gezelter 845 real(kind=dp) :: listSkin = 1.0
422 gezelter 1138
423 mmeineke 377 !! initialize local variables
424 gezelter 1138
425 mmeineke 377 #ifdef IS_MPI
426 chuckv 441 pot_local = 0.0_dp
427 tim 1198 nAtomsInRow = getNatomsInRow(plan_atom_row)
428     nAtomsInCol = getNatomsInCol(plan_atom_col)
429     nGroupsInRow = getNgroupsInRow(plan_group_row)
430     nGroupsInCol = getNgroupsInCol(plan_group_col)
431 mmeineke 377 #else
432     natoms = nlocal
433     #endif
434 gezelter 1138
435 chuckv 900 call doReadyCheck(localError)
436 mmeineke 377 if ( localError .ne. 0 ) then
437 chuckv 900 call handleError("do_force_loop", "Not Initialized")
438 mmeineke 377 error = -1
439     return
440     end if
441     call zero_work_arrays()
442 gezelter 1138
443 mmeineke 377 do_pot = do_pot_c
444     do_stress = do_stress_c
445 gezelter 1138
446 mmeineke 377 ! Gather all information needed by all force loops:
447    
448     #ifdef IS_MPI
449 gezelter 1138
450 tim 1198 call gather(q, q_Row, plan_atom_row_3d)
451     call gather(q, q_Col, plan_atom_col_3d)
452 gezelter 1150
453 tim 1198 call gather(q_group, q_group_Row, plan_group_row_3d)
454     call gather(q_group, q_group_Col, plan_group_col_3d)
455 gezelter 1150
456 chuckv 900 if (FF_UsesDirectionalAtoms() .and. SIM_uses_directional_atoms) then
457 tim 1198 call gather(u_l, u_l_Row, plan_atom_row_3d)
458     call gather(u_l, u_l_Col, plan_atom_col_3d)
459 mmeineke 377
460 tim 1198 call gather(A, A_Row, plan_atom_row_rotation)
461     call gather(A, A_Col, plan_atom_col_rotation)
462 mmeineke 377 endif
463    
464     #endif
465 gezelter 1138
466     !! Begin force loop timing:
467 chuckv 694 #ifdef PROFILE
468     call cpu_time(forceTimeInitial)
469     nloops = nloops + 1
470     #endif
471 gezelter 1138
472 gezelter 1197 loopEnd = PAIR_LOOP
473 chuckv 900 if (FF_RequiresPrepairCalc() .and. SIM_requires_prepair_calc) then
474 gezelter 1197 loopStart = PREPAIR_LOOP
475     else
476     loopStart = PAIR_LOOP
477     endif
478 gezelter 1150
479 gezelter 1197 do loop = loopStart, loopEnd
480 gezelter 1150
481 gezelter 1197 ! See if we need to update neighbor lists
482     ! (but only on the first time through):
483     if (loop .eq. loopStart) then
484 gezelter 1208 #ifdef IS_MPI
485     call checkNeighborList(nGroupsInRow, q_group_row, listSkin, &
486     update_nlist)
487     #else
488     call checkNeighborList(nGroups, q_group, listSkin, &
489     update_nlist)
490     #endif
491 gezelter 1197 endif
492 gezelter 1138
493     if (update_nlist) then
494 gezelter 1197 !! save current configuration and construct neighbor list
495 tim 1198 #ifdef IS_MPI
496 gezelter 1208 call saveNeighborList(nGroupsInRow, q_group_row)
497 tim 1198 #else
498     call saveNeighborList(nGroups, q_group)
499     #endif
500 gezelter 1138 neighborListSize = size(list)
501 gezelter 1197 nlist = 0
502     endif
503    
504     istart = 1
505 gezelter 1192 #ifdef IS_MPI
506 tim 1198 iend = nGroupsInRow
507 gezelter 1192 #else
508 tim 1198 iend = nGroups - 1
509 gezelter 1192 #endif
510 gezelter 1197 outer: do i = istart, iend
511 gezelter 1208
512 gezelter 1197 if (update_nlist) point(i) = nlist + 1
513    
514 tim 1198 n_in_i = groupStartRow(i+1) - groupStartRow(i)
515 gezelter 1197
516     if (update_nlist) then
517 gezelter 1192 #ifdef IS_MPI
518     jstart = 1
519 tim 1198 jend = nGroupsInCol
520 gezelter 1192 #else
521     jstart = i+1
522 tim 1198 jend = nGroups
523 gezelter 1192 #endif
524 gezelter 1197 else
525     jstart = point(i)
526     jend = point(i+1) - 1
527     ! make sure group i has neighbors
528     if (jstart .gt. jend) cycle outer
529     endif
530    
531     do jnab = jstart, jend
532     if (update_nlist) then
533     j = jnab
534     else
535     j = list(jnab)
536     endif
537 gezelter 1199
538 gezelter 1192 #ifdef IS_MPI
539 gezelter 1197 call get_interatomic_vector(q_group_Row(:,i), &
540     q_group_Col(:,j), d_grp, rgrpsq)
541 gezelter 1192 #else
542 gezelter 1197 call get_interatomic_vector(q_group(:,i), &
543     q_group(:,j), d_grp, rgrpsq)
544     #endif
545 gezelter 1199
546 gezelter 1197 if (rgrpsq < rlistsq) then
547     if (update_nlist) then
548 gezelter 1138 nlist = nlist + 1
549    
550     if (nlist > neighborListSize) then
551 tim 1198 #ifdef IS_MPI
552     call expandNeighborList(nGroupsInRow, listerror)
553     #else
554     call expandNeighborList(nGroups, listerror)
555     #endif
556 gezelter 1138 if (listerror /= 0) then
557     error = -1
558     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
559     return
560     end if
561     neighborListSize = size(list)
562     endif
563    
564 gezelter 1197 list(nlist) = j
565     endif
566 gezelter 1138
567 gezelter 1197 if (loop .eq. PAIR_LOOP) then
568     vij = 0.0d0
569     fij(1:3) = 0.0d0
570 mmeineke 377 endif
571    
572 gezelter 1150 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
573     in_switching_region)
574 gezelter 1197
575 tim 1198 n_in_j = groupStartCol(j+1) - groupStartCol(j)
576 gezelter 1150
577 tim 1198 do ia = groupStartRow(i), groupStartRow(i+1)-1
578 gezelter 1197
579 tim 1198 atom1 = groupListRow(ia)
580 gezelter 1150
581 tim 1198 inner: do jb = groupStartCol(j), groupStartCol(j+1)-1
582 gezelter 1150
583 tim 1198 atom2 = groupListCol(jb)
584 gezelter 1197
585     if (skipThisPair(atom1, atom2)) cycle inner
586 gezelter 1214
587 gezelter 1169 if ((n_in_i .eq. 1).and.(n_in_j .eq. 1)) then
588     d_atm(1:3) = d_grp(1:3)
589     ratmsq = rgrpsq
590     else
591 gezelter 1192 #ifdef IS_MPI
592 gezelter 1169 call get_interatomic_vector(q_Row(:,atom1), &
593     q_Col(:,atom2), d_atm, ratmsq)
594 gezelter 1192 #else
595     call get_interatomic_vector(q(:,atom1), &
596     q(:,atom2), d_atm, ratmsq)
597     #endif
598 gezelter 1169 endif
599 gezelter 1214
600 gezelter 1197 if (loop .eq. PREPAIR_LOOP) then
601 gezelter 1192 #ifdef IS_MPI
602 gezelter 1197 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
603     rgrpsq, d_grp, do_pot, do_stress, &
604     u_l, A, f, t, pot_local)
605 gezelter 1192 #else
606 gezelter 1197 call do_prepair(atom1, atom2, ratmsq, d_atm, sw, &
607     rgrpsq, d_grp, do_pot, do_stress, &
608     u_l, A, f, t, pot)
609     #endif
610     else
611     #ifdef IS_MPI
612     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
613     do_pot, &
614     u_l, A, f, t, pot_local, vpair, fpair)
615 gezelter 1192 #else
616 gezelter 1197 call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
617     do_pot, &
618     u_l, A, f, t, pot, vpair, fpair)
619 gezelter 1192 #endif
620 gezelter 1214
621 gezelter 1197 vij = vij + vpair
622     fij(1:3) = fij(1:3) + fpair(1:3)
623     endif
624     enddo inner
625     enddo
626 gezelter 1138
627 gezelter 1197 if (loop .eq. PAIR_LOOP) then
628     if (in_switching_region) then
629     swderiv = vij*dswdr/rgrp
630     fij(1) = fij(1) + swderiv*d_grp(1)
631     fij(2) = fij(2) + swderiv*d_grp(2)
632     fij(3) = fij(3) + swderiv*d_grp(3)
633 gezelter 1150
634 tim 1198 do ia=groupStartRow(i), groupStartRow(i+1)-1
635     atom1=groupListRow(ia)
636     mf = mfactRow(atom1)
637 gezelter 1192 #ifdef IS_MPI
638 gezelter 1197 f_Row(1,atom1) = f_Row(1,atom1) + swderiv*d_grp(1)*mf
639     f_Row(2,atom1) = f_Row(2,atom1) + swderiv*d_grp(2)*mf
640     f_Row(3,atom1) = f_Row(3,atom1) + swderiv*d_grp(3)*mf
641 gezelter 1192 #else
642 gezelter 1197 f(1,atom1) = f(1,atom1) + swderiv*d_grp(1)*mf
643     f(2,atom1) = f(2,atom1) + swderiv*d_grp(2)*mf
644     f(3,atom1) = f(3,atom1) + swderiv*d_grp(3)*mf
645 gezelter 1192 #endif
646 gezelter 1197 enddo
647    
648 tim 1198 do jb=groupStartCol(j), groupStartCol(j+1)-1
649     atom2=groupListCol(jb)
650     mf = mfactCol(atom2)
651 gezelter 1192 #ifdef IS_MPI
652 gezelter 1197 f_Col(1,atom2) = f_Col(1,atom2) - swderiv*d_grp(1)*mf
653     f_Col(2,atom2) = f_Col(2,atom2) - swderiv*d_grp(2)*mf
654     f_Col(3,atom2) = f_Col(3,atom2) - swderiv*d_grp(3)*mf
655 gezelter 1192 #else
656 gezelter 1197 f(1,atom2) = f(1,atom2) - swderiv*d_grp(1)*mf
657     f(2,atom2) = f(2,atom2) - swderiv*d_grp(2)*mf
658     f(3,atom2) = f(3,atom2) - swderiv*d_grp(3)*mf
659 gezelter 1192 #endif
660 gezelter 1197 enddo
661     endif
662 gezelter 1150
663 gezelter 1197 if (do_stress) call add_stress_tensor(d_grp, fij)
664     endif
665     end if
666     enddo
667     enddo outer
668    
669     if (update_nlist) then
670 gezelter 1192 #ifdef IS_MPI
671 tim 1198 point(nGroupsInRow + 1) = nlist + 1
672 gezelter 1197 #else
673 tim 1198 point(nGroups) = nlist + 1
674 gezelter 1192 #endif
675 gezelter 1197 if (loop .eq. PREPAIR_LOOP) then
676     ! we just did the neighbor list update on the first
677     ! pass, so we don't need to do it
678     ! again on the second pass
679     update_nlist = .false.
680 mmeineke 377 endif
681 gezelter 1197 endif
682    
683     if (loop .eq. PREPAIR_LOOP) then
684     call do_preforce(nlocal, pot)
685     endif
686    
687     enddo
688 gezelter 1169
689 gezelter 1138 !! Do timing
690 chuckv 694 #ifdef PROFILE
691     call cpu_time(forceTimeFinal)
692     forceTime = forceTime + forceTimeFinal - forceTimeInitial
693 gezelter 1150 #endif
694 gezelter 1138
695 mmeineke 377 #ifdef IS_MPI
696     !!distribute forces
697 gezelter 1138
698 chuckv 438 f_temp = 0.0_dp
699 tim 1198 call scatter(f_Row,f_temp,plan_atom_row_3d)
700 chuckv 438 do i = 1,nlocal
701     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
702     end do
703 gezelter 1138
704 chuckv 438 f_temp = 0.0_dp
705 tim 1198 call scatter(f_Col,f_temp,plan_atom_col_3d)
706 mmeineke 377 do i = 1,nlocal
707     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
708     end do
709    
710 chuckv 900 if (FF_UsesDirectionalAtoms() .and. SIM_uses_directional_atoms) then
711 chuckv 438 t_temp = 0.0_dp
712 tim 1198 call scatter(t_Row,t_temp,plan_atom_row_3d)
713 chuckv 438 do i = 1,nlocal
714     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
715     end do
716     t_temp = 0.0_dp
717 tim 1198 call scatter(t_Col,t_temp,plan_atom_col_3d)
718 mmeineke 377
719     do i = 1,nlocal
720     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
721     end do
722     endif
723    
724     if (do_pot) then
725     ! scatter/gather pot_row into the members of my column
726 tim 1198 call scatter(pot_Row, pot_Temp, plan_atom_row)
727 gezelter 1138
728 mmeineke 377 ! scatter/gather pot_local into all other procs
729     ! add resultant to get total pot
730     do i = 1, nlocal
731     pot_local = pot_local + pot_Temp(i)
732     enddo
733 chuckv 439
734     pot_Temp = 0.0_DP
735 gezelter 1138
736 tim 1198 call scatter(pot_Col, pot_Temp, plan_atom_col)
737 mmeineke 377 do i = 1, nlocal
738     pot_local = pot_local + pot_Temp(i)
739     enddo
740 gezelter 1150
741 gezelter 1138 endif
742 mmeineke 377 #endif
743 gezelter 1138
744 chuckv 900 if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then
745 mmeineke 377
746 chuckv 900 if (FF_uses_RF .and. SIM_uses_RF) then
747 gezelter 1138
748 mmeineke 377 #ifdef IS_MPI
749 tim 1198 call scatter(rf_Row,rf,plan_atom_row_3d)
750     call scatter(rf_Col,rf_Temp,plan_atom_col_3d)
751 mmeineke 377 do i = 1,nlocal
752     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
753     end do
754     #endif
755    
756 chuckv 898 do i = 1, nLocal
757 gezelter 1138
758 mmeineke 377 rfpot = 0.0_DP
759     #ifdef IS_MPI
760     me_i = atid_row(i)
761     #else
762     me_i = atid(i)
763     #endif
764 gezelter 1138
765 chuckv 900 if (PropertyMap(me_i)%is_DP) then
766 gezelter 1138
767 chuckv 900 mu_i = PropertyMap(me_i)%dipole_moment
768 gezelter 1138
769 mmeineke 377 !! The reaction field needs to include a self contribution
770     !! to the field:
771 chuckv 900 call accumulate_self_rf(i, mu_i, u_l)
772 mmeineke 377 !! Get the reaction field contribution to the
773     !! potential and torques:
774     call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
775     #ifdef IS_MPI
776     pot_local = pot_local + rfpot
777     #else
778     pot = pot + rfpot
779    
780     #endif
781     endif
782     enddo
783     endif
784     endif
785 gezelter 1138
786    
787 mmeineke 377 #ifdef IS_MPI
788 gezelter 1138
789 mmeineke 377 if (do_pot) then
790 chuckv 441 pot = pot + pot_local
791 mmeineke 377 !! we assume the c code will do the allreduce to get the total potential
792     !! we could do it right here if we needed to...
793     endif
794 gezelter 1138
795 mmeineke 377 if (do_stress) then
796 gezelter 1138 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
797 mmeineke 377 mpi_comm_world,mpi_err)
798 chuckv 470 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
799 mmeineke 377 mpi_comm_world,mpi_err)
800     endif
801 gezelter 1138
802 mmeineke 377 #else
803 gezelter 1138
804 mmeineke 377 if (do_stress) then
805     tau = tau_Temp
806     virial = virial_Temp
807     endif
808 mmeineke 887
809 mmeineke 377 #endif
810 gezelter 1197
811 mmeineke 377 end subroutine do_force_loop
812 gezelter 1138
813 gezelter 1192 subroutine do_pair(i, j, rijsq, d, sw, do_pot, &
814     u_l, A, f, t, pot, vpair, fpair)
815 mmeineke 377
816 gezelter 1150 real( kind = dp ) :: pot, vpair, sw
817 gezelter 1192 real( kind = dp ), dimension(3) :: fpair
818 gezelter 1138 real( kind = dp ), dimension(nLocal) :: mfact
819 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
820 gezelter 1138 real( kind = dp ), dimension(9,nLocal) :: A
821     real( kind = dp ), dimension(3,nLocal) :: f
822     real( kind = dp ), dimension(3,nLocal) :: t
823 mmeineke 377
824 gezelter 1192 logical, intent(inout) :: do_pot
825 mmeineke 377 integer, intent(in) :: i, j
826 gezelter 1150 real ( kind = dp ), intent(inout) :: rijsq
827     real ( kind = dp ) :: r
828     real ( kind = dp ), intent(inout) :: d(3)
829 mmeineke 377 integer :: me_i, me_j
830 gezelter 946
831 mmeineke 377 r = sqrt(rijsq)
832 gezelter 1163 vpair = 0.0d0
833 gezelter 1192 fpair(1:3) = 0.0d0
834 mmeineke 377
835     #ifdef IS_MPI
836     me_i = atid_row(i)
837     me_j = atid_col(j)
838     #else
839     me_i = atid(i)
840     me_j = atid(j)
841     #endif
842 chuckv 894
843 chuckv 900 if (FF_uses_LJ .and. SIM_uses_LJ) then
844 gezelter 946
845     if ( PropertyMap(me_i)%is_LJ .and. PropertyMap(me_j)%is_LJ ) then
846 gezelter 1192 call do_lj_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, do_pot)
847 gezelter 946 endif
848    
849 mmeineke 377 endif
850 gezelter 946
851     if (FF_uses_charges .and. SIM_uses_charges) then
852    
853     if (PropertyMap(me_i)%is_Charge .and. PropertyMap(me_j)%is_Charge) then
854 gezelter 1192 call do_charge_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, do_pot)
855 gezelter 946 endif
856    
857     endif
858    
859 chuckv 900 if (FF_uses_dipoles .and. SIM_uses_dipoles) then
860 mmeineke 377
861 chuckv 900 if ( PropertyMap(me_i)%is_DP .and. PropertyMap(me_j)%is_DP) then
862 gezelter 1192 call do_dipole_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, u_l, f, t, &
863     do_pot)
864 chuckv 900 if (FF_uses_RF .and. SIM_uses_RF) then
865 gezelter 1160 call accumulate_rf(i, j, r, u_l, sw)
866 gezelter 1192 call rf_correct_forces(i, j, d, r, u_l, sw, f, fpair)
867 gezelter 946 endif
868 mmeineke 377 endif
869 gezelter 946
870 mmeineke 377 endif
871    
872 chuckv 900 if (FF_uses_Sticky .and. SIM_uses_sticky) then
873 mmeineke 377
874 chuckv 900 if ( PropertyMap(me_i)%is_Sticky .and. PropertyMap(me_j)%is_Sticky) then
875 gezelter 1192 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, A, f, t, &
876     do_pot)
877 mmeineke 377 endif
878 gezelter 946
879 mmeineke 377 endif
880    
881    
882 chuckv 900 if (FF_uses_GB .and. SIM_uses_GB) then
883 mmeineke 377
884 chuckv 900 if ( PropertyMap(me_i)%is_GB .and. PropertyMap(me_j)%is_GB) then
885 gezelter 1192 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, u_l, f, t, &
886     do_pot)
887 mmeineke 377 endif
888 chuckv 894
889 mmeineke 377 endif
890 gezelter 946
891     if (FF_uses_EAM .and. SIM_uses_EAM) then
892    
893     if ( PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) then
894 gezelter 1192 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, &
895     do_pot)
896 gezelter 946 endif
897    
898     endif
899 gezelter 1150
900 mmeineke 377 end subroutine do_pair
901    
902 gezelter 1192 subroutine do_prepair(i, j, rijsq, d, sw, rcijsq, dc, &
903 gezelter 1138 do_pot, do_stress, u_l, A, f, t, pot)
904 gezelter 1192
905     real( kind = dp ) :: pot, sw
906 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
907     real (kind=dp), dimension(9,nLocal) :: A
908     real (kind=dp), dimension(3,nLocal) :: f
909     real (kind=dp), dimension(3,nLocal) :: t
910 chuckv 631
911     logical, intent(inout) :: do_pot, do_stress
912     integer, intent(in) :: i, j
913 gezelter 1138 real ( kind = dp ), intent(inout) :: rijsq, rcijsq
914     real ( kind = dp ) :: r, rc
915     real ( kind = dp ), intent(inout) :: d(3), dc(3)
916 chuckv 631
917     logical :: is_EAM_i, is_EAM_j
918    
919     integer :: me_i, me_j
920    
921 gezelter 1138
922     r = sqrt(rijsq)
923     if (SIM_uses_molecular_cutoffs) then
924     rc = sqrt(rcijsq)
925     else
926     rc = r
927     endif
928 chuckv 631
929 chuckv 669
930 gezelter 1217 #ifdef IS_MPI
931 chuckv 631 me_i = atid_row(i)
932 gezelter 1217 me_j = atid_col(j)
933     #else
934 chuckv 631 me_i = atid(i)
935 gezelter 1217 me_j = atid(j)
936 chuckv 631 #endif
937 gezelter 1192
938 chuckv 900 if (FF_uses_EAM .and. SIM_uses_EAM) then
939 gezelter 1192
940 chuckv 900 if (PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) &
941 chuckv 648 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
942 gezelter 1192
943 chuckv 631 endif
944 chuckv 900
945 chuckv 673 end subroutine do_prepair
946 gezelter 1192
947    
948 gezelter 1138 subroutine do_preforce(nlocal,pot)
949     integer :: nlocal
950     real( kind = dp ) :: pot
951    
952     if (FF_uses_EAM .and. SIM_uses_EAM) then
953     call calc_EAM_preforce_Frho(nlocal,pot)
954     endif
955    
956    
957     end subroutine do_preforce
958    
959    
960     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
961    
962     real (kind = dp), dimension(3) :: q_i
963     real (kind = dp), dimension(3) :: q_j
964     real ( kind = dp ), intent(out) :: r_sq
965     real( kind = dp ) :: d(3), scaled(3)
966     integer i
967    
968     d(1:3) = q_j(1:3) - q_i(1:3)
969    
970     ! Wrap back into periodic box if necessary
971     if ( SIM_uses_PBC ) then
972 mmeineke 377
973 gezelter 1138 if( .not.boxIsOrthorhombic ) then
974     ! calc the scaled coordinates.
975    
976     scaled = matmul(HmatInv, d)
977    
978     ! wrap the scaled coordinates
979    
980     scaled = scaled - anint(scaled)
981    
982    
983     ! calc the wrapped real coordinates from the wrapped scaled
984     ! coordinates
985    
986     d = matmul(Hmat,scaled)
987    
988     else
989     ! calc the scaled coordinates.
990    
991     do i = 1, 3
992     scaled(i) = d(i) * HmatInv(i,i)
993    
994     ! wrap the scaled coordinates
995    
996     scaled(i) = scaled(i) - anint(scaled(i))
997    
998     ! calc the wrapped real coordinates from the wrapped scaled
999     ! coordinates
1000    
1001     d(i) = scaled(i)*Hmat(i,i)
1002     enddo
1003     endif
1004    
1005     endif
1006    
1007     r_sq = dot_product(d,d)
1008    
1009     end subroutine get_interatomic_vector
1010 mmeineke 377
1011 gezelter 1138 subroutine zero_work_arrays()
1012    
1013     #ifdef IS_MPI
1014    
1015     q_Row = 0.0_dp
1016     q_Col = 0.0_dp
1017 mmeineke 377
1018 gezelter 1150 q_group_Row = 0.0_dp
1019     q_group_Col = 0.0_dp
1020 gezelter 1138
1021     u_l_Row = 0.0_dp
1022     u_l_Col = 0.0_dp
1023    
1024     A_Row = 0.0_dp
1025     A_Col = 0.0_dp
1026    
1027     f_Row = 0.0_dp
1028     f_Col = 0.0_dp
1029     f_Temp = 0.0_dp
1030    
1031     t_Row = 0.0_dp
1032     t_Col = 0.0_dp
1033     t_Temp = 0.0_dp
1034    
1035     pot_Row = 0.0_dp
1036     pot_Col = 0.0_dp
1037     pot_Temp = 0.0_dp
1038    
1039     rf_Row = 0.0_dp
1040     rf_Col = 0.0_dp
1041     rf_Temp = 0.0_dp
1042    
1043 mmeineke 377 #endif
1044 chuckv 673
1045 gezelter 1138 if (FF_uses_EAM .and. SIM_uses_EAM) then
1046     call clean_EAM()
1047     endif
1048    
1049     rf = 0.0_dp
1050     tau_Temp = 0.0_dp
1051     virial_Temp = 0.0_dp
1052     end subroutine zero_work_arrays
1053    
1054     function skipThisPair(atom1, atom2) result(skip_it)
1055     integer, intent(in) :: atom1
1056     integer, intent(in), optional :: atom2
1057     logical :: skip_it
1058     integer :: unique_id_1, unique_id_2
1059     integer :: me_i,me_j
1060     integer :: i
1061    
1062     skip_it = .false.
1063    
1064     !! there are a number of reasons to skip a pair or a particle
1065     !! mostly we do this to exclude atoms who are involved in short
1066     !! range interactions (bonds, bends, torsions), but we also need
1067     !! to exclude some overcounted interactions that result from
1068     !! the parallel decomposition
1069    
1070 mmeineke 377 #ifdef IS_MPI
1071 gezelter 1138 !! in MPI, we have to look up the unique IDs for each atom
1072 tim 1198 unique_id_1 = AtomRowToGlobal(atom1)
1073 mmeineke 377 #else
1074 gezelter 1138 !! in the normal loop, the atom numbers are unique
1075     unique_id_1 = atom1
1076 mmeineke 377 #endif
1077 gezelter 1138
1078     !! We were called with only one atom, so just check the global exclude
1079     !! list for this atom
1080     if (.not. present(atom2)) then
1081     do i = 1, nExcludes_global
1082     if (excludesGlobal(i) == unique_id_1) then
1083     skip_it = .true.
1084     return
1085     end if
1086     end do
1087     return
1088     end if
1089    
1090 mmeineke 377 #ifdef IS_MPI
1091 tim 1198 unique_id_2 = AtomColToGlobal(atom2)
1092 mmeineke 377 #else
1093 gezelter 1138 unique_id_2 = atom2
1094 mmeineke 377 #endif
1095 gezelter 1138
1096 mmeineke 377 #ifdef IS_MPI
1097 gezelter 1138 !! this situation should only arise in MPI simulations
1098     if (unique_id_1 == unique_id_2) then
1099     skip_it = .true.
1100     return
1101     end if
1102    
1103     !! this prevents us from doing the pair on multiple processors
1104     if (unique_id_1 < unique_id_2) then
1105     if (mod(unique_id_1 + unique_id_2,2) == 0) then
1106     skip_it = .true.
1107     return
1108     endif
1109     else
1110     if (mod(unique_id_1 + unique_id_2,2) == 1) then
1111     skip_it = .true.
1112     return
1113     endif
1114     endif
1115 mmeineke 377 #endif
1116 gezelter 1138
1117     !! the rest of these situations can happen in all simulations:
1118     do i = 1, nExcludes_global
1119     if ((excludesGlobal(i) == unique_id_1) .or. &
1120     (excludesGlobal(i) == unique_id_2)) then
1121     skip_it = .true.
1122     return
1123     endif
1124     enddo
1125    
1126 tim 1206 do i = 1, nSkipsForAtom(atom1)
1127     if (skipsForAtom(atom1, i) .eq. unique_id_2) then
1128 gezelter 1183 skip_it = .true.
1129     return
1130 gezelter 1138 endif
1131     end do
1132    
1133     return
1134     end function skipThisPair
1135 chuckv 441
1136 gezelter 1138 function FF_UsesDirectionalAtoms() result(doesit)
1137     logical :: doesit
1138     doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
1139     FF_uses_GB .or. FF_uses_RF
1140     end function FF_UsesDirectionalAtoms
1141    
1142     function FF_RequiresPrepairCalc() result(doesit)
1143     logical :: doesit
1144     doesit = FF_uses_EAM
1145     end function FF_RequiresPrepairCalc
1146    
1147     function FF_RequiresPostpairCalc() result(doesit)
1148     logical :: doesit
1149     doesit = FF_uses_RF
1150     end function FF_RequiresPostpairCalc
1151    
1152 chuckv 883 #ifdef PROFILE
1153 gezelter 1138 function getforcetime() result(totalforcetime)
1154     real(kind=dp) :: totalforcetime
1155     totalforcetime = forcetime
1156     end function getforcetime
1157 chuckv 883 #endif
1158 gezelter 1138
1159     !! This cleans componets of force arrays belonging only to fortran
1160 gezelter 1192
1161     subroutine add_stress_tensor(dpair, fpair)
1162    
1163     real( kind = dp ), dimension(3), intent(in) :: dpair, fpair
1164    
1165     ! because the d vector is the rj - ri vector, and
1166     ! because fx, fy, fz are the force on atom i, we need a
1167     ! negative sign here:
1168    
1169     tau_Temp(1) = tau_Temp(1) - dpair(1) * fpair(1)
1170     tau_Temp(2) = tau_Temp(2) - dpair(1) * fpair(2)
1171     tau_Temp(3) = tau_Temp(3) - dpair(1) * fpair(3)
1172     tau_Temp(4) = tau_Temp(4) - dpair(2) * fpair(1)
1173     tau_Temp(5) = tau_Temp(5) - dpair(2) * fpair(2)
1174     tau_Temp(6) = tau_Temp(6) - dpair(2) * fpair(3)
1175     tau_Temp(7) = tau_Temp(7) - dpair(3) * fpair(1)
1176     tau_Temp(8) = tau_Temp(8) - dpair(3) * fpair(2)
1177     tau_Temp(9) = tau_Temp(9) - dpair(3) * fpair(3)
1178    
1179     virial_Temp = virial_Temp + &
1180     (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
1181    
1182     end subroutine add_stress_tensor
1183 gezelter 1138
1184 mmeineke 377 end module do_Forces
1185 gezelter 1192