ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 1163
Committed: Wed May 12 14:30:12 2004 UTC (20 years, 4 months ago) by gezelter
File size: 44122 byte(s)
Log Message:
bug fixes for cutoffGroups

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 1163 !! @version $Id: do_Forces.F90,v 1.56 2004-05-12 14:29:24 gezelter Exp $, $Date: 2004-05-12 14:29:24 $, $Name: not supported by cvs2svn $, $Revision: 1.56 $
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 chuckv 898 real ( kind = dp ), dimension(3,nLocal) :: q
374 gezelter 1138 !! molecular center-of-mass position array
375 gezelter 1150 real ( kind = dp ), dimension(3,nGroup) :: q_group
376 mmeineke 377 !! Rotation Matrix for each long range particle in simulation.
377 chuckv 898 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 1150 real( kind = DP ) :: ratmsq, rgrpsq, rgrp, vpair, vab
406     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 chuckv 897 integer :: me_i, me_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 mmeineke 377 point(i) = nlist + 1
684    
685 gezelter 1150 do j = 1, ncol_group
686 mmeineke 377
687 gezelter 1150 call get_interatomic_vector(q_group_Row(:,i), &
688     q_group_Col(:,j), d_grp, rgrpsq)
689 mmeineke 377
690 gezelter 1150 if (rgrpsq < rlistsq) then
691 mmeineke 377 nlist = nlist + 1
692    
693     if (nlist > neighborListSize) then
694 gezelter 1150 call expandNeighborList(nGroup, listerror)
695 mmeineke 377 if (listerror /= 0) then
696     error = -1
697     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
698     return
699     end if
700     neighborListSize = size(list)
701     endif
702    
703     list(nlist) = j
704 mmeineke 626
705 gezelter 1150 vab = 0.0d0
706     call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
707     in_switching_region)
708    
709     do ia = groupStart(i), groupStart(i+1)-1
710     atom1 = groupList(ia)
711    
712     inner: do jb = groupStart(j), groupStart(j+1)-1
713 gezelter 1163 atom2 = groupList(jb)
714    
715 gezelter 1150 if (skipThisPair(atom1, atom2)) cycle inner
716    
717     call get_interatomic_vector(q_Row(:,atom1), &
718     q_Col(:,atom2), d_atm, ratmsq)
719    
720     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
721     do_pot, do_stress, &
722     u_l, A, f, t, pot_local, vpair)
723    
724     vab = vab + vpair
725     enddo inner
726     enddo
727    
728     if (in_switching_region) then
729     swderiv = vab*dswdr/rgrp
730    
731     do ia=groupStart(i), groupStart(i+1)-1
732     atom1=groupList(ia)
733     mf = mfact(atom1)
734     f_Row(1,atom1) = f_Row(1,atom1) - swderiv*d_grp(1)*mf
735     f_Row(2,atom1) = f_Row(2,atom1) - swderiv*d_grp(2)*mf
736     f_Row(3,atom1) = f_Row(3,atom1) - swderiv*d_grp(3)*mf
737     enddo
738    
739     do jb=groupStart(j), groupStart(j+1)-1
740     atom2=groupList(jb)
741     mf = mfact(atom2)
742     f_Col(1,atom2) = f_Col(1,atom2) + swderiv*d_grp(1)*mf
743     f_Col(2,atom2) = f_Col(2,atom2) + swderiv*d_grp(2)*mf
744     f_Col(3,atom2) = f_Col(3,atom2) + swderiv*d_grp(3)*mf
745     enddo
746 gezelter 1138 endif
747 gezelter 1150
748     end if
749     enddo
750 mmeineke 377 enddo
751 gezelter 1150 point(nrow_group + 1) = nlist + 1
752 gezelter 1138
753 mmeineke 377 else !! (of update_check)
754 gezelter 1138
755 mmeineke 377 ! use the list to find the neighbors
756 gezelter 1150 do i = 1, nrow_group
757 mmeineke 377 JBEG = POINT(i)
758     JEND = POINT(i+1) - 1
759 gezelter 1150 ! check that group i has neighbors
760 mmeineke 377 if (jbeg .le. jend) then
761    
762     do jnab = jbeg, jend
763     j = list(jnab)
764 gezelter 1150
765     call get_interatomic_vector(q_group_Row(:,i), &
766     q_group_Col(:,j), d_grp, rgrpsq)
767 gezelter 1138
768 gezelter 1150 vab = 0.0d0
769     call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
770     in_switching_region)
771 gezelter 1138
772 gezelter 1150 do ia = groupStart(i), groupStart(i+1)-1
773     atom1 = groupList(ia)
774    
775     do jb = groupStart(j), groupStart(j+1)-1
776     atom2 = groupList(jb)
777 gezelter 1163 write(*,*) 'doing ', atom1, atom2
778 gezelter 1150
779     call get_interatomic_vector(q_Row(:,atom1), &
780     q_Col(:,atom2), d_atm, ratmsq)
781    
782     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
783     do_pot, do_stress, &
784     u_l, A, f, t, pot_local, vpair)
785    
786     vab = vab + vpair
787    
788     enddo
789     enddo
790 gezelter 1138
791 gezelter 1150 if (in_switching_region) then
792     swderiv = vab*dswdr/rgrp
793    
794     do ia=groupStart(i), groupStart(i+1)-1
795     atom1=groupList(ia)
796     mf = mfact(atom1)
797     f_Row(1,atom1) = f_Row(1,atom1) - swderiv*d_grp(1)*mf
798     f_Row(2,atom1) = f_Row(2,atom1) - swderiv*d_grp(2)*mf
799     f_Row(3,atom1) = f_Row(3,atom1) - swderiv*d_grp(3)*mf
800     enddo
801    
802     do jb=groupStart(j), groupStart(j+1)-1
803     atom2=groupList(jb)
804     mf = mfact(atom2)
805     f_Col(1,atom2) = f_Col(1,atom2) + swderiv*d_grp(1)*mf
806     f_Col(2,atom2) = f_Col(2,atom2) + swderiv*d_grp(2)*mf
807     f_Col(3,atom2) = f_Col(3,atom2) + swderiv*d_grp(3)*mf
808     enddo
809     endif
810    
811 mmeineke 377 enddo
812     endif
813     enddo
814     endif
815    
816     #else
817 gezelter 1150
818 mmeineke 377 if (update_nlist) then
819 gezelter 1138
820 gezelter 1150 !! save current configuration, construct neighbor list,
821     !! and calculate forces
822 mmeineke 377
823 gezelter 1150 call saveNeighborList(nGroup, q_group)
824    
825 mmeineke 377 neighborListSize = size(list)
826 gezelter 1150 nlist = 0
827 gezelter 1138
828 gezelter 1150 do i = 1, nGroup-1
829 mmeineke 377 point(i) = nlist + 1
830    
831 gezelter 1150 do j = i+1, nGroup
832 mmeineke 377
833 gezelter 1150 call get_interatomic_vector(q_group(:,i), &
834     q_group(:,j), d_grp, rgrpsq)
835 gezelter 1138
836 gezelter 1150 if (rgrpsq < rlistsq) then
837 mmeineke 377 nlist = nlist + 1
838 gezelter 1138
839 mmeineke 377 if (nlist > neighborListSize) then
840 gezelter 1150 call expandNeighborList(nGroup, listerror)
841 mmeineke 377 if (listerror /= 0) then
842     error = -1
843     write(DEFAULT_ERROR,*) "ERROR: nlist > list size and max allocations exceeded."
844     return
845     end if
846     neighborListSize = size(list)
847     endif
848    
849     list(nlist) = j
850 gezelter 1150
851     vab = 0.0d0
852     call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
853     in_switching_region)
854 mmeineke 377
855 gezelter 1150 do ia = groupStart(i), groupStart(i+1)-1
856     atom1 = groupList(ia)
857    
858     inner: do jb = groupStart(j), groupStart(j+1)-1
859     atom2 = groupList(jb)
860    
861     if (skipThisPair(atom1, atom2)) cycle inner
862    
863     call get_interatomic_vector(q(:,atom1), &
864     q(:,atom2), d_atm, ratmsq)
865    
866     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
867     do_pot, do_stress, &
868     u_l, A, f, t, pot, vpair)
869    
870     vab = vab + vpair
871    
872     enddo inner
873     enddo
874    
875     if (in_switching_region) then
876     swderiv = vab*dswdr/rgrp
877     do ia=groupStart(i), groupStart(i+1)-1
878     atom1=groupList(ia)
879     mf = mfact(atom1)
880     f(1,atom1) = f(1,atom1) - swderiv*d_grp(1)*mf
881     f(2,atom1) = f(2,atom1) - swderiv*d_grp(2)*mf
882     f(3,atom1) = f(3,atom1) - swderiv*d_grp(3)*mf
883     enddo
884    
885     do jb=groupStart(j), groupStart(j+1)-1
886     atom2=groupList(jb)
887     mf = mfact(atom2)
888     f(1,atom2) = f(1,atom2) + swderiv*d_grp(1)*mf
889     f(2,atom2) = f(2,atom2) + swderiv*d_grp(2)*mf
890     f(3,atom2) = f(3,atom2) + swderiv*d_grp(3)*mf
891     enddo
892 gezelter 1138 endif
893 gezelter 1150
894     end if
895     enddo
896 mmeineke 377 enddo
897 gezelter 1150 point(nGroup) = nlist + 1
898 mmeineke 377
899 gezelter 1150 else !! (of update_check)
900 mmeineke 377
901     ! use the list to find the neighbors
902 gezelter 1150 do i = 1, nGroup-1
903 mmeineke 377 JBEG = POINT(i)
904     JEND = POINT(i+1) - 1
905 gezelter 1150 ! check that group i has neighbors
906 mmeineke 377 if (jbeg .le. jend) then
907    
908     do jnab = jbeg, jend
909     j = list(jnab)
910    
911 gezelter 1150 call get_interatomic_vector(q_group(:,i), &
912     q_group(:,j), d_grp, rgrpsq)
913 gezelter 1138
914 gezelter 1163 vab = 0.0d0
915    
916 gezelter 1150 call get_switch(rgrpsq, sw, dswdr, rgrp, group_switch, &
917     in_switching_region)
918    
919     do ia = groupStart(i), groupStart(i+1)-1
920     atom1 = groupList(ia)
921    
922     do jb = groupStart(j), groupStart(j+1)-1
923     atom2 = groupList(jb)
924    
925     call get_interatomic_vector(q(:,atom1), &
926     q(:,atom2), d_atm, ratmsq)
927    
928     call do_pair(atom1, atom2, ratmsq, d_atm, sw, &
929     do_pot, do_stress, &
930     u_l, A, f, t, pot, vpair)
931    
932     vab = vab + vpair
933 gezelter 1163
934 gezelter 1150 enddo
935     enddo
936    
937     if (in_switching_region) then
938     swderiv = vab*dswdr/rgrp
939    
940     do ia=groupStart(i), groupStart(i+1)-1
941     atom1=groupList(ia)
942     mf = mfact(atom1)
943     f(1,atom1) = f(1,atom1) - swderiv*d_grp(1)*mf
944     f(2,atom1) = f(2,atom1) - swderiv*d_grp(2)*mf
945     f(3,atom1) = f(3,atom1) - swderiv*d_grp(3)*mf
946     enddo
947    
948     do jb=groupStart(j), groupStart(j+1)-1
949     atom2=groupList(jb)
950     mf = mfact(atom2)
951     f(1,atom2) = f(1,atom2) + swderiv*d_grp(1)*mf
952     f(2,atom2) = f(2,atom2) + swderiv*d_grp(2)*mf
953     f(3,atom2) = f(3,atom2) + swderiv*d_grp(3)*mf
954     enddo
955 gezelter 1138 endif
956 mmeineke 377 enddo
957     endif
958     enddo
959     endif
960 gezelter 1150
961 mmeineke 377 #endif
962    
963     ! phew, done with main loop.
964 gezelter 1138
965     !! Do timing
966 chuckv 694 #ifdef PROFILE
967     call cpu_time(forceTimeFinal)
968     forceTime = forceTime + forceTimeFinal - forceTimeInitial
969 gezelter 1150 #endif
970 gezelter 1138
971 mmeineke 377 #ifdef IS_MPI
972     !!distribute forces
973 gezelter 1138
974 chuckv 438 f_temp = 0.0_dp
975     call scatter(f_Row,f_temp,plan_row3d)
976     do i = 1,nlocal
977     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
978     end do
979 gezelter 1138
980 chuckv 438 f_temp = 0.0_dp
981 mmeineke 377 call scatter(f_Col,f_temp,plan_col3d)
982     do i = 1,nlocal
983     f(1:3,i) = f(1:3,i) + f_temp(1:3,i)
984     end do
985    
986 chuckv 900 if (FF_UsesDirectionalAtoms() .and. SIM_uses_directional_atoms) then
987 chuckv 438 t_temp = 0.0_dp
988     call scatter(t_Row,t_temp,plan_row3d)
989     do i = 1,nlocal
990     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
991     end do
992     t_temp = 0.0_dp
993 mmeineke 377 call scatter(t_Col,t_temp,plan_col3d)
994    
995     do i = 1,nlocal
996     t(1:3,i) = t(1:3,i) + t_temp(1:3,i)
997     end do
998     endif
999    
1000     if (do_pot) then
1001     ! scatter/gather pot_row into the members of my column
1002     call scatter(pot_Row, pot_Temp, plan_row)
1003 gezelter 1138
1004 mmeineke 377 ! scatter/gather pot_local into all other procs
1005     ! add resultant to get total pot
1006     do i = 1, nlocal
1007     pot_local = pot_local + pot_Temp(i)
1008     enddo
1009 chuckv 439
1010     pot_Temp = 0.0_DP
1011 gezelter 1138
1012 mmeineke 377 call scatter(pot_Col, pot_Temp, plan_col)
1013     do i = 1, nlocal
1014     pot_local = pot_local + pot_Temp(i)
1015     enddo
1016 gezelter 1150
1017 gezelter 1138 endif
1018 mmeineke 377 #endif
1019 gezelter 1138
1020 chuckv 900 if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then
1021 mmeineke 377
1022 chuckv 900 if (FF_uses_RF .and. SIM_uses_RF) then
1023 gezelter 1138
1024 mmeineke 377 #ifdef IS_MPI
1025     call scatter(rf_Row,rf,plan_row3d)
1026     call scatter(rf_Col,rf_Temp,plan_col3d)
1027     do i = 1,nlocal
1028     rf(1:3,i) = rf(1:3,i) + rf_Temp(1:3,i)
1029     end do
1030     #endif
1031    
1032 chuckv 898 do i = 1, nLocal
1033 gezelter 1138
1034 mmeineke 377 rfpot = 0.0_DP
1035     #ifdef IS_MPI
1036     me_i = atid_row(i)
1037     #else
1038     me_i = atid(i)
1039     #endif
1040 gezelter 1138
1041 chuckv 900 if (PropertyMap(me_i)%is_DP) then
1042 gezelter 1138
1043 chuckv 900 mu_i = PropertyMap(me_i)%dipole_moment
1044 gezelter 1138
1045 mmeineke 377 !! The reaction field needs to include a self contribution
1046     !! to the field:
1047 chuckv 900 call accumulate_self_rf(i, mu_i, u_l)
1048 mmeineke 377 !! Get the reaction field contribution to the
1049     !! potential and torques:
1050     call reaction_field_final(i, mu_i, u_l, rfpot, t, do_pot)
1051     #ifdef IS_MPI
1052     pot_local = pot_local + rfpot
1053     #else
1054     pot = pot + rfpot
1055    
1056     #endif
1057     endif
1058     enddo
1059     endif
1060     endif
1061 gezelter 1138
1062    
1063 mmeineke 377 #ifdef IS_MPI
1064 gezelter 1138
1065 mmeineke 377 if (do_pot) then
1066 chuckv 441 pot = pot + pot_local
1067 mmeineke 377 !! we assume the c code will do the allreduce to get the total potential
1068     !! we could do it right here if we needed to...
1069     endif
1070 gezelter 1138
1071 mmeineke 377 if (do_stress) then
1072 gezelter 1138 call mpi_allreduce(tau_Temp, tau, 9,mpi_double_precision,mpi_sum, &
1073 mmeineke 377 mpi_comm_world,mpi_err)
1074 chuckv 470 call mpi_allreduce(virial_Temp, virial,1,mpi_double_precision,mpi_sum, &
1075 mmeineke 377 mpi_comm_world,mpi_err)
1076     endif
1077 gezelter 1138
1078 mmeineke 377 #else
1079 gezelter 1138
1080 mmeineke 377 if (do_stress) then
1081     tau = tau_Temp
1082     virial = virial_Temp
1083     endif
1084 mmeineke 887
1085 mmeineke 377 #endif
1086 mmeineke 887
1087    
1088 mmeineke 377 end subroutine do_force_loop
1089 gezelter 1150
1090 gezelter 1138
1091 gezelter 1150 subroutine do_pair(i, j, rijsq, d, sw, do_pot, do_stress, &
1092     u_l, A, f, t, pot, vpair)
1093 mmeineke 377
1094 gezelter 1150 real( kind = dp ) :: pot, vpair, sw
1095 gezelter 1138 real( kind = dp ), dimension(nLocal) :: mfact
1096 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
1097 gezelter 1138 real( kind = dp ), dimension(9,nLocal) :: A
1098     real( kind = dp ), dimension(3,nLocal) :: f
1099     real( kind = dp ), dimension(3,nLocal) :: t
1100 mmeineke 377
1101     logical, intent(inout) :: do_pot, do_stress
1102     integer, intent(in) :: i, j
1103 gezelter 1150 real ( kind = dp ), intent(inout) :: rijsq
1104     real ( kind = dp ) :: r
1105     real ( kind = dp ), intent(inout) :: d(3)
1106 mmeineke 377 integer :: me_i, me_j
1107 gezelter 946
1108 mmeineke 377 r = sqrt(rijsq)
1109 gezelter 1163 vpair = 0.0d0
1110 mmeineke 377
1111     #ifdef IS_MPI
1112 gezelter 490 if (tagRow(i) .eq. tagColumn(j)) then
1113     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
1114     endif
1115 mmeineke 377 me_i = atid_row(i)
1116     me_j = atid_col(j)
1117     #else
1118     me_i = atid(i)
1119     me_j = atid(j)
1120     #endif
1121 chuckv 894
1122 chuckv 900 if (FF_uses_LJ .and. SIM_uses_LJ) then
1123 gezelter 946
1124     if ( PropertyMap(me_i)%is_LJ .and. PropertyMap(me_j)%is_LJ ) then
1125 gezelter 1163 !write(*,*) 'calling lj with'
1126     !write(*,*) i, j, r, rijsq
1127     !write(*,'(3es12.3)') d(1), d(2), d(3)
1128     !write(*,'(3es12.3)') sw, vpair, pot
1129     !write(*,*)
1130    
1131 gezelter 1150 call do_lj_pair(i, j, d, r, rijsq, sw, vpair, pot, f, do_pot, &
1132     do_stress)
1133 gezelter 946 endif
1134    
1135 mmeineke 377 endif
1136 gezelter 946
1137     if (FF_uses_charges .and. SIM_uses_charges) then
1138    
1139     if (PropertyMap(me_i)%is_Charge .and. PropertyMap(me_j)%is_Charge) then
1140 gezelter 1150 call do_charge_pair(i, j, d, r, rijsq, sw, vpair, pot, f, do_pot, &
1141     do_stress)
1142 gezelter 946 endif
1143    
1144     endif
1145    
1146 chuckv 900 if (FF_uses_dipoles .and. SIM_uses_dipoles) then
1147 mmeineke 377
1148 chuckv 900 if ( PropertyMap(me_i)%is_DP .and. PropertyMap(me_j)%is_DP) then
1149 gezelter 1150 call do_dipole_pair(i, j, d, r, rijsq, sw, vpair, pot, u_l, f, t, &
1150 mmeineke 377 do_pot, do_stress)
1151 chuckv 900 if (FF_uses_RF .and. SIM_uses_RF) then
1152 gezelter 1160 call accumulate_rf(i, j, r, u_l, sw)
1153     call rf_correct_forces(i, j, d, r, u_l, sw, f, do_stress)
1154 gezelter 946 endif
1155 mmeineke 377 endif
1156 gezelter 946
1157 mmeineke 377 endif
1158    
1159 chuckv 900 if (FF_uses_Sticky .and. SIM_uses_sticky) then
1160 mmeineke 377
1161 chuckv 900 if ( PropertyMap(me_i)%is_Sticky .and. PropertyMap(me_j)%is_Sticky) then
1162 gezelter 1150 call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, pot, A, f, t, &
1163 mmeineke 377 do_pot, do_stress)
1164     endif
1165 gezelter 946
1166 mmeineke 377 endif
1167    
1168    
1169 chuckv 900 if (FF_uses_GB .and. SIM_uses_GB) then
1170 mmeineke 377
1171 chuckv 900 if ( PropertyMap(me_i)%is_GB .and. PropertyMap(me_j)%is_GB) then
1172 gezelter 1150 call do_gb_pair(i, j, d, r, rijsq, sw, vpair, pot, u_l, f, t, &
1173 mmeineke 377 do_pot, do_stress)
1174     endif
1175 chuckv 894
1176 mmeineke 377 endif
1177 gezelter 946
1178     if (FF_uses_EAM .and. SIM_uses_EAM) then
1179    
1180     if ( PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) then
1181 gezelter 1150 call do_eam_pair(i, j, d, r, rijsq, sw, vpair, pot, f, &
1182     do_pot, do_stress)
1183 gezelter 946 endif
1184    
1185     endif
1186 gezelter 1150
1187 mmeineke 377 end subroutine do_pair
1188    
1189 gezelter 1138 subroutine do_prepair(i, j, rijsq, d, rcijsq, dc, &
1190     do_pot, do_stress, u_l, A, f, t, pot)
1191 chuckv 631 real( kind = dp ) :: pot
1192 chuckv 898 real( kind = dp ), dimension(3,nLocal) :: u_l
1193     real (kind=dp), dimension(9,nLocal) :: A
1194     real (kind=dp), dimension(3,nLocal) :: f
1195     real (kind=dp), dimension(3,nLocal) :: t
1196 chuckv 631
1197     logical, intent(inout) :: do_pot, do_stress
1198     integer, intent(in) :: i, j
1199 gezelter 1138 real ( kind = dp ), intent(inout) :: rijsq, rcijsq
1200     real ( kind = dp ) :: r, rc
1201     real ( kind = dp ), intent(inout) :: d(3), dc(3)
1202 chuckv 631
1203     logical :: is_EAM_i, is_EAM_j
1204    
1205     integer :: me_i, me_j
1206    
1207 gezelter 1138
1208     r = sqrt(rijsq)
1209     if (SIM_uses_molecular_cutoffs) then
1210     rc = sqrt(rcijsq)
1211     else
1212     rc = r
1213     endif
1214 chuckv 631
1215 chuckv 669
1216 chuckv 631 #ifdef IS_MPI
1217     if (tagRow(i) .eq. tagColumn(j)) then
1218     write(0,*) 'do_pair is doing', i , j, tagRow(i), tagColumn(j)
1219     endif
1220    
1221     me_i = atid_row(i)
1222     me_j = atid_col(j)
1223    
1224     #else
1225    
1226     me_i = atid(i)
1227     me_j = atid(j)
1228    
1229     #endif
1230 chuckv 673
1231 chuckv 900 if (FF_uses_EAM .and. SIM_uses_EAM) then
1232    
1233     if (PropertyMap(me_i)%is_EAM .and. PropertyMap(me_j)%is_EAM) &
1234 chuckv 648 call calc_EAM_prepair_rho(i, j, d, r, rijsq )
1235 chuckv 900
1236 chuckv 631 endif
1237 chuckv 900
1238 chuckv 673 end subroutine do_prepair
1239 chuckv 648
1240    
1241    
1242 chuckv 673
1243 gezelter 1138 subroutine do_preforce(nlocal,pot)
1244     integer :: nlocal
1245     real( kind = dp ) :: pot
1246    
1247     if (FF_uses_EAM .and. SIM_uses_EAM) then
1248     call calc_EAM_preforce_Frho(nlocal,pot)
1249     endif
1250    
1251    
1252     end subroutine do_preforce
1253    
1254    
1255     subroutine get_interatomic_vector(q_i, q_j, d, r_sq)
1256    
1257     real (kind = dp), dimension(3) :: q_i
1258     real (kind = dp), dimension(3) :: q_j
1259     real ( kind = dp ), intent(out) :: r_sq
1260     real( kind = dp ) :: d(3), scaled(3)
1261     integer i
1262    
1263     d(1:3) = q_j(1:3) - q_i(1:3)
1264    
1265     ! Wrap back into periodic box if necessary
1266     if ( SIM_uses_PBC ) then
1267 mmeineke 377
1268 gezelter 1138 if( .not.boxIsOrthorhombic ) then
1269     ! calc the scaled coordinates.
1270    
1271     scaled = matmul(HmatInv, d)
1272    
1273     ! wrap the scaled coordinates
1274    
1275     scaled = scaled - anint(scaled)
1276    
1277    
1278     ! calc the wrapped real coordinates from the wrapped scaled
1279     ! coordinates
1280    
1281     d = matmul(Hmat,scaled)
1282    
1283     else
1284     ! calc the scaled coordinates.
1285    
1286     do i = 1, 3
1287     scaled(i) = d(i) * HmatInv(i,i)
1288    
1289     ! wrap the scaled coordinates
1290    
1291     scaled(i) = scaled(i) - anint(scaled(i))
1292    
1293     ! calc the wrapped real coordinates from the wrapped scaled
1294     ! coordinates
1295    
1296     d(i) = scaled(i)*Hmat(i,i)
1297     enddo
1298     endif
1299    
1300     endif
1301    
1302     r_sq = dot_product(d,d)
1303    
1304     end subroutine get_interatomic_vector
1305 mmeineke 377
1306 gezelter 1138 subroutine zero_work_arrays()
1307    
1308     #ifdef IS_MPI
1309    
1310     q_Row = 0.0_dp
1311     q_Col = 0.0_dp
1312 mmeineke 377
1313 gezelter 1150 q_group_Row = 0.0_dp
1314     q_group_Col = 0.0_dp
1315 gezelter 1138
1316     u_l_Row = 0.0_dp
1317     u_l_Col = 0.0_dp
1318    
1319     A_Row = 0.0_dp
1320     A_Col = 0.0_dp
1321    
1322     f_Row = 0.0_dp
1323     f_Col = 0.0_dp
1324     f_Temp = 0.0_dp
1325    
1326     t_Row = 0.0_dp
1327     t_Col = 0.0_dp
1328     t_Temp = 0.0_dp
1329    
1330     pot_Row = 0.0_dp
1331     pot_Col = 0.0_dp
1332     pot_Temp = 0.0_dp
1333    
1334     rf_Row = 0.0_dp
1335     rf_Col = 0.0_dp
1336     rf_Temp = 0.0_dp
1337    
1338 mmeineke 377 #endif
1339 chuckv 673
1340 gezelter 1138 if (FF_uses_EAM .and. SIM_uses_EAM) then
1341     call clean_EAM()
1342     endif
1343    
1344     rf = 0.0_dp
1345     tau_Temp = 0.0_dp
1346     virial_Temp = 0.0_dp
1347     end subroutine zero_work_arrays
1348    
1349     function skipThisPair(atom1, atom2) result(skip_it)
1350     integer, intent(in) :: atom1
1351     integer, intent(in), optional :: atom2
1352     logical :: skip_it
1353     integer :: unique_id_1, unique_id_2
1354     integer :: me_i,me_j
1355     integer :: i
1356    
1357     skip_it = .false.
1358    
1359     !! there are a number of reasons to skip a pair or a particle
1360     !! mostly we do this to exclude atoms who are involved in short
1361     !! range interactions (bonds, bends, torsions), but we also need
1362     !! to exclude some overcounted interactions that result from
1363     !! the parallel decomposition
1364    
1365 mmeineke 377 #ifdef IS_MPI
1366 gezelter 1138 !! in MPI, we have to look up the unique IDs for each atom
1367     unique_id_1 = tagRow(atom1)
1368 mmeineke 377 #else
1369 gezelter 1138 !! in the normal loop, the atom numbers are unique
1370     unique_id_1 = atom1
1371 mmeineke 377 #endif
1372 gezelter 1138
1373     !! We were called with only one atom, so just check the global exclude
1374     !! list for this atom
1375     if (.not. present(atom2)) then
1376     do i = 1, nExcludes_global
1377     if (excludesGlobal(i) == unique_id_1) then
1378     skip_it = .true.
1379     return
1380     end if
1381     end do
1382     return
1383     end if
1384    
1385 mmeineke 377 #ifdef IS_MPI
1386 gezelter 1138 unique_id_2 = tagColumn(atom2)
1387 mmeineke 377 #else
1388 gezelter 1138 unique_id_2 = atom2
1389 mmeineke 377 #endif
1390 gezelter 1138
1391 mmeineke 377 #ifdef IS_MPI
1392 gezelter 1138 !! this situation should only arise in MPI simulations
1393     if (unique_id_1 == unique_id_2) then
1394     skip_it = .true.
1395     return
1396     end if
1397    
1398     !! this prevents us from doing the pair on multiple processors
1399     if (unique_id_1 < unique_id_2) then
1400     if (mod(unique_id_1 + unique_id_2,2) == 0) then
1401     skip_it = .true.
1402     return
1403     endif
1404     else
1405     if (mod(unique_id_1 + unique_id_2,2) == 1) then
1406     skip_it = .true.
1407     return
1408     endif
1409     endif
1410 mmeineke 377 #endif
1411 gezelter 1138
1412     !! the rest of these situations can happen in all simulations:
1413     do i = 1, nExcludes_global
1414     if ((excludesGlobal(i) == unique_id_1) .or. &
1415     (excludesGlobal(i) == unique_id_2)) then
1416     skip_it = .true.
1417     return
1418     endif
1419     enddo
1420    
1421     do i = 1, nExcludes_local
1422     if (excludesLocal(1,i) == unique_id_1) then
1423     if (excludesLocal(2,i) == unique_id_2) then
1424     skip_it = .true.
1425     return
1426     endif
1427     else
1428     if (excludesLocal(1,i) == unique_id_2) then
1429     if (excludesLocal(2,i) == unique_id_1) then
1430     skip_it = .true.
1431     return
1432     endif
1433     endif
1434     endif
1435     end do
1436    
1437     return
1438     end function skipThisPair
1439 chuckv 441
1440 gezelter 1138 function FF_UsesDirectionalAtoms() result(doesit)
1441     logical :: doesit
1442     doesit = FF_uses_dipoles .or. FF_uses_sticky .or. &
1443     FF_uses_GB .or. FF_uses_RF
1444     end function FF_UsesDirectionalAtoms
1445    
1446     function FF_RequiresPrepairCalc() result(doesit)
1447     logical :: doesit
1448     doesit = FF_uses_EAM
1449     end function FF_RequiresPrepairCalc
1450    
1451     function FF_RequiresPostpairCalc() result(doesit)
1452     logical :: doesit
1453     doesit = FF_uses_RF
1454     end function FF_RequiresPostpairCalc
1455    
1456 chuckv 883 #ifdef PROFILE
1457 gezelter 1138 function getforcetime() result(totalforcetime)
1458     real(kind=dp) :: totalforcetime
1459     totalforcetime = forcetime
1460     end function getforcetime
1461 chuckv 883 #endif
1462 gezelter 1138
1463     !! This cleans componets of force arrays belonging only to fortran
1464    
1465 mmeineke 377 end module do_Forces