--- trunk/OOPSE-4/src/UseTheForce/doForces.F90 2005/08/09 22:33:37 2270 +++ trunk/OOPSE-4/src/UseTheForce/doForces.F90 2005/09/18 20:45:38 2309 @@ -45,7 +45,7 @@ !! @author Charles F. Vardeman II !! @author Matthew Meineke -!! @version $Id: doForces.F90,v 1.28 2005-08-09 22:33:37 gezelter Exp $, $Date: 2005-08-09 22:33:37 $, $Name: not supported by cvs2svn $, $Revision: 1.28 $ +!! @version $Id: doForces.F90,v 1.46 2005-09-18 20:45:38 chrisfen Exp $, $Date: 2005-09-18 20:45:38 $, $Name: not supported by cvs2svn $, $Revision: 1.46 $ 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,8 +73,10 @@ 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 @@ -83,24 +85,31 @@ module doForces 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_Dipoles logical, save :: FF_uses_GayBerne logical, save :: FF_uses_EAM - logical, save :: FF_uses_RF logical, save :: SIM_uses_DirectionalAtoms logical, save :: SIM_uses_EAM - logical, save :: SIM_uses_RF logical, save :: SIM_requires_postpair_calc logical, save :: SIM_requires_prepair_calc logical, save :: SIM_uses_PBC + integer, save :: electrostaticSummationMethod + public :: init_FF + public :: setDefaultCutoffs public :: do_force_loop public :: createInteractionHash public :: createGtypeCutoffMap + public :: getStickyCut + public :: getStickyPowerCut + public :: getGayBerneCut + public :: getEAMCut + public :: getShapeCut #ifdef PROFILE public :: getforcetime @@ -122,6 +131,9 @@ module doForces real(kind=dp) :: rlistsq end type gtypeCutoffs type(gtypeCutoffs), dimension(:,:), allocatable :: gtypeCutoffMap + + integer, save :: cutoffPolicy = TRADITIONAL_CUTOFF_POLICY + real(kind=dp),save :: defaultRcut, defaultRsw, defaultRlist contains @@ -146,7 +158,8 @@ contains logical :: j_is_GB logical :: j_is_EAM logical :: j_is_Shape - + real(kind=dp) :: myRcut + status = 0 if (.not. associated(atypes)) then @@ -179,13 +192,6 @@ contains call getElementProperty(atypes, i, "is_EAM", i_is_EAM) call getElementProperty(atypes, i, "is_Shape", i_is_Shape) - if (i_is_LJ) then - thisCut = getDefaultLJCutoff(i) - if (thisCut .gt. atypeMaxCutoff(i)) atypeMaxCutoff(i) = thisCut - endif - - - do j = i, nAtypes iHash = 0 @@ -238,13 +244,24 @@ contains haveInteractionHash = .true. end subroutine createInteractionHash - subroutine createGtypeCutoffMap(defaultRcut, defaultSkinThickness, stat) - - real(kind=dp), intent(in), optional :: defaultRCut, defaultSkinThickness - integer, intent(out) :: stat + subroutine createGtypeCutoffMap(stat) - integer :: myStatus, nAtypes + 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) @@ -254,24 +271,167 @@ contains 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 - - atypeMaxCutoff(i) = - + 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 - haveGtypeCutoffMap = .true. + 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 + 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_EAM = SimUsesEAM() - SIM_uses_RF = SimUsesRF() SIM_requires_postpair_calc = SimRequiresPostpairCalc() SIM_requires_prepair_calc = SimRequiresPrepairCalc() SIM_uses_PBC = SimUsesPBC() @@ -312,11 +472,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!' @@ -341,10 +501,10 @@ contains end subroutine doReadyCheck - subroutine init_FF(use_RF_c, thisStat) - - logical, intent(in) :: use_RF_c + subroutine init_FF(thisESM, thisStat) + integer, intent(in) :: thisESM + real(kind=dp), intent(in) :: dampingAlpha integer, intent(out) :: thisStat integer :: my_status, nMatches integer, pointer :: MatchList(:) => null() @@ -353,8 +513,7 @@ contains !! assume things are copacetic, unless they aren't thisStat = 0 - !! Fortran's version of a cast: - FF_uses_RF = use_RF_c + electrostaticSummationMethod = thisESM !! init_FF is called *after* all of the atom types have been !! defined in atype_module using the new_atype subroutine. @@ -385,15 +544,15 @@ contains haveSaneForceField = .true. - !! check to make sure the FF_uses_RF setting makes sense + !! check to make sure the reaction field setting makes sense if (FF_uses_Dipoles) then - if (FF_uses_RF) then + if (electrostaticSummationMethod == 3) then dielect = getDielect() call initialize_rf(dielect) endif else - if (FF_uses_RF) then + if (electrostaticSummationMethod == 3) then write(default_error,*) 'Using Reaction Field with no dipoles? Huh?' thisStat = -1 haveSaneForceField = .false. @@ -575,12 +734,6 @@ contains #endif outer: do i = istart, iend -#ifdef IS_MPI - me_i = atid_row(i) -#else - me_i = atid(i) -#endif - if (update_nlist) point(i) = nlist + 1 n_in_i = groupStartRow(i+1) - groupStartRow(i) @@ -617,7 +770,7 @@ contains q_group(:,j), d_grp, rgrpsq) #endif - if (rgrpsq < InteractionHash(me_i,me_j)%rListsq) then + if (rgrpsq < gtypeCutoffMap(groupToGtype(i),groupToGtype(j))%rListsq) then if (update_nlist) then nlist = nlist + 1 @@ -817,7 +970,7 @@ contains if (FF_RequiresPostpairCalc() .and. SIM_requires_postpair_calc) then - if (FF_uses_RF .and. SIM_uses_RF) then + if (electrostaticSummationMethod == 3) then #ifdef IS_MPI call scatter(rf_Row,rf,plan_atom_row_3d) @@ -901,7 +1054,6 @@ 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 :: iHash @@ -928,7 +1080,7 @@ contains call doElectrostaticPair(i, j, d, r, rijsq, sw, vpair, fpair, & pot, eFrame, f, t, do_pot) - if (FF_uses_RF .and. SIM_uses_RF) then + if (electrostaticSummationMethod == 3) then ! CHECK ME (RF needs to know about all electrostatic types) call accumulate_rf(i, j, r, eFrame, sw) @@ -991,6 +1143,8 @@ contains integer :: me_i, me_j, iHash + r = sqrt(rijsq) + #ifdef IS_MPI me_i = atid_row(i) me_j = atid_col(j) @@ -1208,7 +1362,7 @@ contains function FF_RequiresPostpairCalc() result(doesit) logical :: doesit - doesit = FF_uses_RF + if (electrostaticSummationMethod == 3) doesit = .true. end function FF_RequiresPostpairCalc #ifdef PROFILE