ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 900
Committed: Tue Jan 6 18:54:57 2004 UTC (20 years, 6 months ago) by chuckv
File size: 32217 byte(s)
Log Message:
Making do_Forces a little more sane

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