ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/do_Forces.F90
Revision: 946
Committed: Thu Jan 15 04:33:24 2004 UTC (20 years, 5 months ago) by gezelter
File size: 32948 byte(s)
Log Message:
changes for charge charge interactions

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