ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/doForces.F90
Revision: 1628
Committed: Thu Oct 21 20:15:31 2004 UTC (19 years, 8 months ago) by gezelter
File size: 34256 byte(s)
Log Message:
Breaky Breaky.   Fixey Fixey.

File Contents

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