--- trunk/OOPSE-4/src/UseTheForce/doForces.F90 2005/06/27 21:01:36 2259 +++ trunk/OOPSE-4/src/UseTheForce/doForces.F90 2005/09/15 22:05:21 2301 @@ -45,7 +45,7 @@ !! @author Charles F. Vardeman II !! @author Matthew Meineke -!! @version $Id: doForces.F90,v 1.20 2005-06-27 21:01:30 gezelter Exp $, $Date: 2005-06-27 21:01:30 $, $Name: not supported by cvs2svn $, $Revision: 1.20 $ +!! @version $Id: doForces.F90,v 1.44 2005-09-15 22:05:17 gezelter Exp $, $Date: 2005-09-15 22:05:17 $, $Name: not supported by cvs2svn $, $Revision: 1.44 $ module doForces @@ -58,7 +58,7 @@ module doForces use lj use sticky use electrostatic_module - use reaction_field + use reaction_field_module use gb_pair use shapes use vector_class @@ -73,54 +73,45 @@ module doForces #define __FORTRAN90 #include "UseTheForce/fSwitchingFunction.h" +#include "UseTheForce/fCutoffPolicy.h" #include "UseTheForce/DarkSide/fInteractionMap.h" + INTEGER, PARAMETER:: PREPAIR_LOOP = 1 INTEGER, PARAMETER:: PAIR_LOOP = 2 - logical, save :: haveRlist = .false. logical, save :: haveNeighborList = .false. logical, save :: haveSIMvariables = .false. - logical, save :: havePropertyMap = .false. logical, save :: haveSaneForceField = .false. + logical, save :: haveInteractionHash = .false. + logical, save :: haveGtypeCutoffMap = .false. + logical, save :: haveRlist = .false. logical, save :: FF_uses_DirectionalAtoms - logical, save :: FF_uses_LennardJones - logical, save :: FF_uses_Electrostatics - logical, save :: FF_uses_Charges logical, save :: FF_uses_Dipoles - logical, save :: FF_uses_Quadrupoles - logical, save :: FF_uses_Sticky - logical, save :: FF_uses_StickyPower logical, save :: FF_uses_GayBerne logical, save :: FF_uses_EAM - logical, save :: FF_uses_Shapes - logical, save :: FF_uses_FLARB logical, save :: FF_uses_RF logical, save :: SIM_uses_DirectionalAtoms - logical, save :: SIM_uses_LennardJones - logical, save :: SIM_uses_Electrostatics - logical, save :: SIM_uses_Charges - logical, save :: SIM_uses_Dipoles - logical, save :: SIM_uses_Quadrupoles - logical, save :: SIM_uses_Sticky - logical, save :: SIM_uses_StickyPower - logical, save :: SIM_uses_GayBerne logical, save :: SIM_uses_EAM - logical, save :: SIM_uses_Shapes - logical, save :: SIM_uses_FLARB logical, save :: SIM_uses_RF logical, save :: SIM_requires_postpair_calc logical, save :: SIM_requires_prepair_calc logical, save :: SIM_uses_PBC - logical, save :: SIM_uses_molecular_cutoffs - real(kind=dp), save :: rlist, rlistsq + integer, save :: corrMethod public :: init_FF + public :: setDefaultCutoffs public :: do_force_loop - public :: setRlistDF + public :: createInteractionHash + public :: createGtypeCutoffMap + public :: getStickyCut + public :: getStickyPowerCut + public :: getGayBerneCut + public :: getEAMCut + public :: getShapeCut #ifdef PROFILE public :: getforcetime @@ -128,115 +119,327 @@ module doForces real :: forceTimeInitial, forceTimeFinal integer :: nLoops #endif + + !! Variables for cutoff mapping and interaction mapping + ! Bit hash to determine pair-pair interactions. + integer, dimension(:,:), allocatable :: InteractionHash + real(kind=dp), dimension(:), allocatable :: atypeMaxCutoff + real(kind=dp), dimension(:), allocatable :: groupMaxCutoff + integer, dimension(:), allocatable :: groupToGtype + real(kind=dp), dimension(:), allocatable :: gtypeMaxCutoff + type ::gtypeCutoffs + real(kind=dp) :: rcut + real(kind=dp) :: rcutsq + real(kind=dp) :: rlistsq + end type gtypeCutoffs + type(gtypeCutoffs), dimension(:,:), allocatable :: gtypeCutoffMap - type :: Properties - logical :: is_Directional = .false. - logical :: is_LennardJones = .false. - logical :: is_Electrostatic = .false. - logical :: is_Charge = .false. - logical :: is_Dipole = .false. - logical :: is_Quadrupole = .false. - logical :: is_Sticky = .false. - logical :: is_StickyPower = .false. - logical :: is_GayBerne = .false. - logical :: is_EAM = .false. - logical :: is_Shape = .false. - logical :: is_FLARB = .false. - end type Properties - - type(Properties), dimension(:),allocatable :: PropertyMap - + integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY + real(kind=dp),save :: defaultRcut, defaultRsw, defaultRlist + real(kind=dp),save :: rcuti + contains - subroutine setRlistDF( this_rlist ) - - real(kind=dp) :: this_rlist - - rlist = this_rlist - rlistsq = rlist * rlist - - haveRlist = .true. - - end subroutine setRlistDF - - subroutine createPropertyMap(status) + subroutine createInteractionHash(status) integer :: nAtypes - integer :: status + integer, intent(out) :: status integer :: i - logical :: thisProperty - real (kind=DP) :: thisDPproperty + integer :: j + integer :: iHash + !! Test Types + logical :: i_is_LJ + logical :: i_is_Elect + logical :: i_is_Sticky + logical :: i_is_StickyP + logical :: i_is_GB + logical :: i_is_EAM + logical :: i_is_Shape + logical :: j_is_LJ + logical :: j_is_Elect + logical :: j_is_Sticky + logical :: j_is_StickyP + logical :: j_is_GB + logical :: j_is_EAM + logical :: j_is_Shape + real(kind=dp) :: myRcut - status = 0 + status = 0 + if (.not. associated(atypes)) then + call handleError("atype", "atypes was not present before call of createInteractionHash!") + status = -1 + return + endif + nAtypes = getSize(atypes) - + if (nAtypes == 0) then status = -1 return end if - if (.not. allocated(PropertyMap)) then - allocate(PropertyMap(nAtypes)) + if (.not. allocated(InteractionHash)) then + allocate(InteractionHash(nAtypes,nAtypes)) endif + if (.not. allocated(atypeMaxCutoff)) then + allocate(atypeMaxCutoff(nAtypes)) + endif + do i = 1, nAtypes - call getElementProperty(atypes, i, "is_Directional", thisProperty) - PropertyMap(i)%is_Directional = thisProperty + call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ) + call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect) + call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky) + call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP) + call getElementProperty(atypes, i, "is_GayBerne", i_is_GB) + call getElementProperty(atypes, i, "is_EAM", i_is_EAM) + call getElementProperty(atypes, i, "is_Shape", i_is_Shape) - call getElementProperty(atypes, i, "is_LennardJones", thisProperty) - PropertyMap(i)%is_LennardJones = thisProperty + do j = i, nAtypes - call getElementProperty(atypes, i, "is_Electrostatic", thisProperty) - PropertyMap(i)%is_Electrostatic = thisProperty + iHash = 0 + myRcut = 0.0_dp - call getElementProperty(atypes, i, "is_Charge", thisProperty) - PropertyMap(i)%is_Charge = thisProperty + call getElementProperty(atypes, j, "is_LennardJones", j_is_LJ) + call getElementProperty(atypes, j, "is_Electrostatic", j_is_Elect) + call getElementProperty(atypes, j, "is_Sticky", j_is_Sticky) + call getElementProperty(atypes, j, "is_StickyPower", j_is_StickyP) + call getElementProperty(atypes, j, "is_GayBerne", j_is_GB) + call getElementProperty(atypes, j, "is_EAM", j_is_EAM) + call getElementProperty(atypes, j, "is_Shape", j_is_Shape) - call getElementProperty(atypes, i, "is_Dipole", thisProperty) - PropertyMap(i)%is_Dipole = thisProperty + if (i_is_LJ .and. j_is_LJ) then + iHash = ior(iHash, LJ_PAIR) + endif + + if (i_is_Elect .and. j_is_Elect) then + iHash = ior(iHash, ELECTROSTATIC_PAIR) + endif + + if (i_is_Sticky .and. j_is_Sticky) then + iHash = ior(iHash, STICKY_PAIR) + endif - call getElementProperty(atypes, i, "is_Quadrupole", thisProperty) - PropertyMap(i)%is_Quadrupole = thisProperty + if (i_is_StickyP .and. j_is_StickyP) then + iHash = ior(iHash, STICKYPOWER_PAIR) + endif - call getElementProperty(atypes, i, "is_Sticky", thisProperty) - PropertyMap(i)%is_Sticky = thisProperty - - call getElementProperty(atypes, i, "is_StickyPower", thisProperty) - PropertyMap(i)%is_StickyPower = thisProperty + if (i_is_EAM .and. j_is_EAM) then + iHash = ior(iHash, EAM_PAIR) + endif - call getElementProperty(atypes, i, "is_GayBerne", thisProperty) - PropertyMap(i)%is_GayBerne = thisProperty + if (i_is_GB .and. j_is_GB) iHash = ior(iHash, GAYBERNE_PAIR) + if (i_is_GB .and. j_is_LJ) iHash = ior(iHash, GAYBERNE_LJ) + if (i_is_LJ .and. j_is_GB) iHash = ior(iHash, GAYBERNE_LJ) - call getElementProperty(atypes, i, "is_EAM", thisProperty) - PropertyMap(i)%is_EAM = thisProperty + if (i_is_Shape .and. j_is_Shape) iHash = ior(iHash, SHAPE_PAIR) + if (i_is_Shape .and. j_is_LJ) iHash = ior(iHash, SHAPE_LJ) + if (i_is_LJ .and. j_is_Shape) iHash = ior(iHash, SHAPE_LJ) - call getElementProperty(atypes, i, "is_Shape", thisProperty) - PropertyMap(i)%is_Shape = thisProperty - call getElementProperty(atypes, i, "is_FLARB", thisProperty) - PropertyMap(i)%is_FLARB = thisProperty + InteractionHash(i,j) = iHash + InteractionHash(j,i) = iHash + + end do + end do - havePropertyMap = .true. + haveInteractionHash = .true. + end subroutine createInteractionHash - end subroutine createPropertyMap + subroutine createGtypeCutoffMap(stat) + + integer, intent(out), optional :: stat + logical :: i_is_LJ + logical :: i_is_Elect + logical :: i_is_Sticky + logical :: i_is_StickyP + logical :: i_is_GB + logical :: i_is_EAM + logical :: i_is_Shape + logical :: GtypeFound + + integer :: myStatus, nAtypes, i, j, istart, iend, jstart, jend + integer :: n_in_i, me_i, ia, g, atom1, nGroupTypes + integer :: nGroupsInRow + real(kind=dp):: thisSigma, bigSigma, thisRcut, tol, skin + real(kind=dp) :: biggestAtypeCutoff + + stat = 0 + if (.not. haveInteractionHash) then + call createInteractionHash(myStatus) + if (myStatus .ne. 0) then + write(default_error, *) 'createInteractionHash failed in doForces!' + stat = -1 + return + endif + endif +#ifdef IS_MPI + nGroupsInRow = getNgroupsInRow(plan_group_row) +#endif + nAtypes = getSize(atypes) +! Set all of the initial cutoffs to zero. + atypeMaxCutoff = 0.0_dp + do i = 1, nAtypes + if (SimHasAtype(i)) then + call getElementProperty(atypes, i, "is_LennardJones", i_is_LJ) + call getElementProperty(atypes, i, "is_Electrostatic", i_is_Elect) + call getElementProperty(atypes, i, "is_Sticky", i_is_Sticky) + call getElementProperty(atypes, i, "is_StickyPower", i_is_StickyP) + call getElementProperty(atypes, i, "is_GayBerne", i_is_GB) + call getElementProperty(atypes, i, "is_EAM", i_is_EAM) + call getElementProperty(atypes, i, "is_Shape", i_is_Shape) + + + if (i_is_LJ) then + thisRcut = getSigma(i) * 2.5_dp + if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut + endif + if (i_is_Elect) then + thisRcut = defaultRcut + if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut + endif + if (i_is_Sticky) then + thisRcut = getStickyCut(i) + if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut + endif + if (i_is_StickyP) then + thisRcut = getStickyPowerCut(i) + if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut + endif + if (i_is_GB) then + thisRcut = getGayBerneCut(i) + if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut + endif + if (i_is_EAM) then + thisRcut = getEAMCut(i) + if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut + endif + if (i_is_Shape) then + thisRcut = getShapeCut(i) + if (thisRCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisRCut + endif + + if (atypeMaxCutoff(i).gt.biggestAtypeCutoff) then + biggestAtypeCutoff = atypeMaxCutoff(i) + endif + endif + enddo + + nGroupTypes = 0 + + istart = 1 +#ifdef IS_MPI + iend = nGroupsInRow +#else + iend = nGroups +#endif + + !! allocate the groupToGtype and gtypeMaxCutoff here. + if(.not.allocated(groupToGtype)) then + allocate(groupToGtype(iend)) + allocate(groupMaxCutoff(iend)) + allocate(gtypeMaxCutoff(iend)) + groupMaxCutoff = 0.0_dp + gtypeMaxCutoff = 0.0_dp + endif + !! first we do a single loop over the cutoff groups to find the + !! largest cutoff for any atypes present in this group. We also + !! create gtypes at this point. + + tol = 1.0d-6 + + do i = istart, iend + n_in_i = groupStartRow(i+1) - groupStartRow(i) + groupMaxCutoff(i) = 0.0_dp + do ia = groupStartRow(i), groupStartRow(i+1)-1 + atom1 = groupListRow(ia) +#ifdef IS_MPI + me_i = atid_row(atom1) +#else + me_i = atid(atom1) +#endif + if (atypeMaxCutoff(me_i).gt.groupMaxCutoff(i)) then + groupMaxCutoff(i)=atypeMaxCutoff(me_i) + endif + enddo + if (nGroupTypes.eq.0) then + nGroupTypes = nGroupTypes + 1 + gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i) + groupToGtype(i) = nGroupTypes + else + GtypeFound = .false. + do g = 1, nGroupTypes + if ( abs(groupMaxCutoff(i) - gtypeMaxCutoff(g)).lt.tol) then + groupToGtype(i) = g + GtypeFound = .true. + endif + enddo + if (.not.GtypeFound) then + nGroupTypes = nGroupTypes + 1 + gtypeMaxCutoff(nGroupTypes) = groupMaxCutoff(i) + groupToGtype(i) = nGroupTypes + endif + endif + enddo + + !! allocate the gtypeCutoffMap here. + allocate(gtypeCutoffMap(nGroupTypes,nGroupTypes)) + !! then we do a double loop over all the group TYPES to find the cutoff + !! map between groups of two types + + do i = 1, nGroupTypes + do j = 1, nGroupTypes + + select case(cutoffPolicy) + case(TRADITIONAL_CUTOFF_POLICY) + thisRcut = maxval(gtypeMaxCutoff) + case(MIX_CUTOFF_POLICY) + thisRcut = 0.5_dp * (gtypeMaxCutoff(i) + gtypeMaxCutoff(j)) + case(MAX_CUTOFF_POLICY) + thisRcut = max(gtypeMaxCutoff(i), gtypeMaxCutoff(j)) + case default + call handleError("createGtypeCutoffMap", "Unknown Cutoff Policy") + return + end select + gtypeCutoffMap(i,j)%rcut = thisRcut + gtypeCutoffMap(i,j)%rcutsq = thisRcut*thisRcut + skin = defaultRlist - defaultRcut + gtypeCutoffMap(i,j)%rlistsq = (thisRcut + skin)**2 + + enddo + enddo + + haveGtypeCutoffMap = .true. + end subroutine createGtypeCutoffMap + + subroutine setDefaultCutoffs(defRcut, defRsw, defRlist, cutPolicy) + real(kind=dp),intent(in) :: defRcut, defRsw, defRlist + integer, intent(in) :: cutPolicy + + defaultRcut = defRcut + defaultRsw = defRsw + defaultRlist = defRlist + cutoffPolicy = cutPolicy + rcuti = 1.0_dp / defaultRcut + end subroutine setDefaultCutoffs + + subroutine setCutoffPolicy(cutPolicy) + + integer, intent(in) :: cutPolicy + cutoffPolicy = cutPolicy + call createGtypeCutoffMap() + end subroutine setCutoffPolicy + + subroutine setSimVariables() SIM_uses_DirectionalAtoms = SimUsesDirectionalAtoms() - SIM_uses_LennardJones = SimUsesLennardJones() - SIM_uses_Electrostatics = SimUsesElectrostatics() - SIM_uses_Charges = SimUsesCharges() - SIM_uses_Dipoles = SimUsesDipoles() - SIM_uses_Sticky = SimUsesSticky() - SIM_uses_StickyPower = SimUsesStickyPower() - SIM_uses_GayBerne = SimUsesGayBerne() SIM_uses_EAM = SimUsesEAM() - SIM_uses_Shapes = SimUsesShapes() - SIM_uses_FLARB = SimUsesFLARB() - SIM_uses_RF = SimUsesRF() SIM_requires_postpair_calc = SimRequiresPostpairCalc() SIM_requires_prepair_calc = SimRequiresPrepairCalc() SIM_uses_PBC = SimUsesPBC() + SIM_uses_RF = SimUsesRF() haveSIMvariables = .true. @@ -250,14 +453,21 @@ contains error = 0 - if (.not. havePropertyMap) then + if (.not. haveInteractionHash) then + myStatus = 0 + call createInteractionHash(myStatus) + if (myStatus .ne. 0) then + write(default_error, *) 'createInteractionHash failed in doForces!' + error = -1 + return + endif + endif - myStatus = 0 - - call createPropertyMap(myStatus) - + if (.not. haveGtypeCutoffMap) then + myStatus = 0 + call createGtypeCutoffMap(myStatus) if (myStatus .ne. 0) then - write(default_error, *) 'createPropertyMap failed in doForces!' + write(default_error, *) 'createGtypeCutoffMap failed in doForces!' error = -1 return endif @@ -267,11 +477,11 @@ contains call setSimVariables() endif - if (.not. haveRlist) then - write(default_error, *) 'rList has not been set in doForces!' - error = -1 - return - endif + ! if (.not. haveRlist) then + ! write(default_error, *) 'rList has not been set in doForces!' + ! error = -1 + ! return + ! endif if (.not. haveNeighborList) then write(default_error, *) 'neighbor list has not been initialized in doForces!' @@ -296,10 +506,11 @@ contains end subroutine doReadyCheck - subroutine init_FF(use_RF_c, thisStat) + subroutine init_FF(use_RF, correctionMethod, dampingAlpha, thisStat) - logical, intent(in) :: use_RF_c - + logical, intent(in) :: use_RF + integer, intent(in) :: correctionMethod + real(kind=dp), intent(in) :: dampingAlpha integer, intent(out) :: thisStat integer :: my_status, nMatches integer, pointer :: MatchList(:) => null() @@ -309,8 +520,9 @@ contains thisStat = 0 !! Fortran's version of a cast: - FF_uses_RF = use_RF_c + FF_uses_RF = use_RF + !! init_FF is called *after* all of the atom types have been !! defined in atype_module using the new_atype subroutine. !! @@ -318,101 +530,37 @@ contains !! interactions are used by the force field. FF_uses_DirectionalAtoms = .false. - FF_uses_LennardJones = .false. - FF_uses_Electrostatics = .false. - FF_uses_Charges = .false. FF_uses_Dipoles = .false. - FF_uses_Sticky = .false. - FF_uses_StickyPower = .false. FF_uses_GayBerne = .false. FF_uses_EAM = .false. - FF_uses_Shapes = .false. - FF_uses_FLARB = .false. call getMatchingElementList(atypes, "is_Directional", .true., & nMatches, MatchList) if (nMatches .gt. 0) FF_uses_DirectionalAtoms = .true. - call getMatchingElementList(atypes, "is_LennardJones", .true., & - nMatches, MatchList) - if (nMatches .gt. 0) FF_uses_LennardJones = .true. - - call getMatchingElementList(atypes, "is_Electrostatic", .true., & - nMatches, MatchList) - if (nMatches .gt. 0) then - FF_uses_Electrostatics = .true. - endif - - call getMatchingElementList(atypes, "is_Charge", .true., & - nMatches, MatchList) - if (nMatches .gt. 0) then - FF_uses_Charges = .true. - FF_uses_Electrostatics = .true. - endif - call getMatchingElementList(atypes, "is_Dipole", .true., & nMatches, MatchList) - if (nMatches .gt. 0) then - FF_uses_Dipoles = .true. - FF_uses_Electrostatics = .true. - FF_uses_DirectionalAtoms = .true. - endif - - call getMatchingElementList(atypes, "is_Quadrupole", .true., & - nMatches, MatchList) - if (nMatches .gt. 0) then - FF_uses_Quadrupoles = .true. - FF_uses_Electrostatics = .true. - FF_uses_DirectionalAtoms = .true. - endif - - call getMatchingElementList(atypes, "is_Sticky", .true., nMatches, & - MatchList) - if (nMatches .gt. 0) then - FF_uses_Sticky = .true. - FF_uses_DirectionalAtoms = .true. - endif - - call getMatchingElementList(atypes, "is_StickyPower", .true., nMatches, & - MatchList) - if (nMatches .gt. 0) then - FF_uses_StickyPower = .true. - FF_uses_DirectionalAtoms = .true. - endif + if (nMatches .gt. 0) FF_uses_Dipoles = .true. call getMatchingElementList(atypes, "is_GayBerne", .true., & nMatches, MatchList) - if (nMatches .gt. 0) then - FF_uses_GayBerne = .true. - FF_uses_DirectionalAtoms = .true. - endif + if (nMatches .gt. 0) FF_uses_GayBerne = .true. call getMatchingElementList(atypes, "is_EAM", .true., nMatches, MatchList) if (nMatches .gt. 0) FF_uses_EAM = .true. - call getMatchingElementList(atypes, "is_Shape", .true., & - nMatches, MatchList) - if (nMatches .gt. 0) then - FF_uses_Shapes = .true. - FF_uses_DirectionalAtoms = .true. - endif - call getMatchingElementList(atypes, "is_FLARB", .true., & - nMatches, MatchList) - if (nMatches .gt. 0) FF_uses_FLARB = .true. - - !! Assume sanity (for the sake of argument) haveSaneForceField = .true. !! check to make sure the FF_uses_RF setting makes sense - if (FF_uses_dipoles) then + if (FF_uses_Dipoles) then if (FF_uses_RF) then dielect = getDielect() call initialize_rf(dielect) endif else - if (FF_uses_RF) then + if ((corrMethod == 3) .or. FF_uses_RF) then write(default_error,*) 'Using Reaction Field with no dipoles? Huh?' thisStat = -1 haveSaneForceField = .false. @@ -420,16 +568,6 @@ contains endif endif - !sticky module does not contain check_sticky_FF anymore - !if (FF_uses_sticky) then - ! call check_sticky_FF(my_status) - ! if (my_status /= 0) then - ! thisStat = -1 - ! haveSaneForceField = .false. - ! return - ! end if - !endif - if (FF_uses_EAM) then call init_EAM_FF(my_status) if (my_status /= 0) then @@ -449,9 +587,6 @@ contains endif endif - if (FF_uses_GayBerne .and. FF_uses_LennardJones) then - endif - if (.not. haveNeighborList) then !! Create neighbor lists call expandNeighborList(nLocal, my_status) @@ -515,7 +650,7 @@ contains integer :: localError integer :: propPack_i, propPack_j integer :: loopStart, loopEnd, loop - + integer :: iHash real(kind=dp) :: listSkin = 1.0 !! initialize local variables @@ -634,14 +769,16 @@ contains endif #ifdef IS_MPI + me_j = atid_col(j) call get_interatomic_vector(q_group_Row(:,i), & q_group_Col(:,j), d_grp, rgrpsq) #else + me_j = atid(j) call get_interatomic_vector(q_group(:,i), & q_group(:,j), d_grp, rgrpsq) #endif - if (rgrpsq < rlistsq) then + if (rgrpsq < gtypeCutoffMap(groupToGtype(i),groupToGtype(j))%rListsq) then if (update_nlist) then nlist = nlist + 1 @@ -841,7 +978,7 @@ contains if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then - if (FF_uses_RF .and. SIM_uses_RF) then + if ((FF_uses_RF .and. SIM_uses_RF) .or. (corrMethod == 3)) then #ifdef IS_MPI call scatter(rf_Row,rf,plan_atom_row_3d) @@ -859,8 +996,9 @@ contains #else me_i = atid(i) #endif - - if (PropertyMap(me_i)%is_Dipole) then + iHash = InteractionHash(me_i,me_j) + + if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then mu_i = getDipoleMoment(me_i) @@ -924,10 +1062,9 @@ contains real ( kind = dp ), intent(inout) :: rijsq real ( kind = dp ) :: r real ( kind = dp ), intent(inout) :: d(3) - real ( kind = dp ) :: ebalance integer :: me_i, me_j - integer :: iMap + integer :: iHash r = sqrt(rijsq) vpair = 0.0d0 @@ -941,58 +1078,56 @@ contains me_j = atid(j) #endif - iMap = InteractionMap(me_i, me_j)%InteractionHash + iHash = InteractionHash(me_i, me_j) - if ( iand(iMap, LJ_PAIR).ne.0 ) then + if ( iand(iHash, LJ_PAIR).ne.0 ) then call do_lj_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, do_pot) endif - if ( iand(iMap, ELECTROSTATIC_PAIR).ne.0 ) then + if ( iand(iHash, ELECTROSTATIC_PAIR).ne.0 ) then call doElectrostaticPair(i, j, d, r, rijsq, sw, vpair, fpair, & - pot, eFrame, f, t, do_pot) + pot, eFrame, f, t, do_pot, corrMethod, rcuti) - if (FF_uses_dipoles .and. SIM_uses_dipoles) then - if ( PropertyMap(me_i)%is_Dipole .and. & - PropertyMap(me_j)%is_Dipole) then - if (FF_uses_RF .and. SIM_uses_RF) then - call accumulate_rf(i, j, r, eFrame, sw) - call rf_correct_forces(i, j, d, r, eFrame, sw, f, fpair) - endif - endif + if ((FF_uses_RF .and. SIM_uses_RF) .or. (corrMethod == 3)) then + + ! CHECK ME (RF needs to know about all electrostatic types) + call accumulate_rf(i, j, r, eFrame, sw) + call rf_correct_forces(i, j, d, r, eFrame, sw, f, fpair) endif + endif - if ( iand(iMap, STICKY_PAIR).ne.0 ) then + if ( iand(iHash, STICKY_PAIR).ne.0 ) then call do_sticky_pair(i, j, d, r, rijsq, sw, vpair, fpair, & pot, A, f, t, do_pot) endif - if ( iand(iMap, STICKYPOWER_PAIR).ne.0 ) then + if ( iand(iHash, STICKYPOWER_PAIR).ne.0 ) then call do_sticky_power_pair(i, j, d, r, rijsq, sw, vpair, fpair, & pot, A, f, t, do_pot) endif - if ( iand(iMap, GAYBERNE_PAIR).ne.0 ) then + if ( iand(iHash, GAYBERNE_PAIR).ne.0 ) then call do_gb_pair(i, j, d, r, rijsq, sw, vpair, fpair, & pot, A, f, t, do_pot) endif - if ( iand(iMap, GAYBERNE_LJ).ne.0 ) then - call do_gblj_pair(i, j, d, r, rijsq, sw, vpair, fpair, & - pot, A, f, t, do_pot) + if ( iand(iHash, GAYBERNE_LJ).ne.0 ) then + ! call do_gblj_pair(i, j, d, r, rijsq, sw, vpair, fpair, & + ! pot, A, f, t, do_pot) endif - if ( iand(iMap, EAM_PAIR).ne.0 ) then + if ( iand(iHash, EAM_PAIR).ne.0 ) then call do_eam_pair(i, j, d, r, rijsq, sw, vpair, fpair, pot, f, & do_pot) endif - if ( iand(iMap, SHAPE_PAIR).ne.0 ) then + if ( iand(iHash, SHAPE_PAIR).ne.0 ) then call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, & pot, A, f, t, do_pot) endif - if ( iand(iMap, SHAPE_LJ).ne.0 ) then + if ( iand(iHash, SHAPE_LJ).ne.0 ) then call do_shape_pair(i, j, d, r, rijsq, sw, vpair, fpair, & pot, A, f, t, do_pot) endif @@ -1014,8 +1149,10 @@ contains real ( kind = dp ) :: r, rc real ( kind = dp ), intent(inout) :: d(3), dc(3) - integer :: me_i, me_j, iMap + integer :: me_i, me_j, iHash + r = sqrt(rijsq) + #ifdef IS_MPI me_i = atid_row(i) me_j = atid_col(j) @@ -1024,9 +1161,9 @@ contains me_j = atid(j) #endif - iMap = InteractionMap(me_i, me_j)%InteractionHash + iHash = InteractionHash(me_i, me_j) - if ( iand(iMap, EAM_PAIR).ne.0 ) then + if ( iand(iHash, EAM_PAIR).ne.0 ) then call calc_EAM_prepair_rho(i, j, d, r, rijsq ) endif @@ -1223,9 +1360,7 @@ contains function FF_UsesDirectionalAtoms() result(doesit) logical :: doesit - doesit = FF_uses_DirectionalAtoms .or. FF_uses_Dipoles .or. & - FF_uses_Quadrupoles .or. FF_uses_Sticky .or. & - FF_uses_StickyPower .or. FF_uses_GayBerne .or. FF_uses_Shapes + doesit = FF_uses_DirectionalAtoms end function FF_UsesDirectionalAtoms function FF_RequiresPrepairCalc() result(doesit) @@ -1236,6 +1371,7 @@ contains function FF_RequiresPostpairCalc() result(doesit) logical :: doesit doesit = FF_uses_RF + if (corrMethod == 3) doesit = .true. end function FF_RequiresPostpairCalc #ifdef PROFILE