--- trunk/OOPSE/libmdtools/simulation_module.F90 2003/04/11 18:46:37 491 +++ trunk/OOPSE/libmdtools/simulation_module.F90 2004/05/21 15:58:48 1183 @@ -6,6 +6,7 @@ module simulation use force_globals use vector_class use atype_module + use switcheroo #ifdef IS_MPI use mpiSimulation #endif @@ -15,37 +16,35 @@ module simulation #define __FORTRAN90 #include "fSimulation.h" +#include "fSwitchingFunction.h" - type (simtype), public :: thisSim + type (simtype), public, save :: thisSim logical, save :: simulation_setup_complete = .false. integer, public, save :: nLocal, nGlobal + integer, public, save :: nGroup, nGroupGlobal integer, public, save :: nExcludes_Global = 0 integer, public, save :: nExcludes_Local = 0 integer, allocatable, dimension(:,:), public :: excludesLocal - integer, allocatable, dimension(:), public :: excludesGlobal - integer, allocatable, dimension(:), public :: molMembershipList + integer, allocatable, dimension(:), public :: excludesGlobal + integer, allocatable, dimension(:), public :: molMembershipList + integer, allocatable, dimension(:), public :: groupList + integer, allocatable, dimension(:), public :: groupStart + integer, allocatable, dimension(:), public :: nSkipsForAtom + integer, allocatable, dimension(:,:), public :: skipsForAtom + real(kind=dp), allocatable, dimension(:), public :: mfact - real(kind=dp), save :: rcut2 = 0.0_DP - real(kind=dp), save :: rcut6 = 0.0_DP - real(kind=dp), save :: rlist2 = 0.0_DP - real(kind=dp), public, dimension(3), save :: box - + real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv + logical, public, save :: boxIsOrthorhombic public :: SimulationSetup public :: getNlocal public :: setBox - public :: setBox_3d - public :: getBox - public :: setRcut - public :: getRcut - public :: getRlist - public :: getRrf - public :: getRt public :: getDielect public :: SimUsesPBC public :: SimUsesLJ + public :: SimUsesCharges public :: SimUsesDipoles public :: SimUsesSticky public :: SimUsesRF @@ -54,22 +53,12 @@ module simulation public :: SimRequiresPrepairCalc public :: SimRequiresPostpairCalc public :: SimUsesDirectionalAtoms - - interface getBox - module procedure getBox_3d - module procedure getBox_1d - end interface - interface setBox - module procedure setBox_3d - module procedure setBox_1d - end interface - contains subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, & CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, & - CmolMembership, & + CmolMembership, Cmfact, CnGroup, CgroupList, CgroupStart, & status) type (simtype) :: setThisSim @@ -83,7 +72,16 @@ contains integer, dimension(CnGlobal),intent(in) :: CmolMembership !! Result status, success = 0, status = -1 integer, intent(out) :: status - integer :: i, me, thisStat, alloc_stat, myNode + integer :: i, j, me, thisStat, alloc_stat, myNode, id1, id2 + integer :: gStart, gEnd, groupSize, biggestGroupSize, ia + + !! mass factors used for molecular cutoffs + real ( kind = dp ), dimension(CnLocal) :: Cmfact + integer, intent(in):: CnGroup + integer, dimension(CnLocal),intent(in) :: CgroupList + integer, dimension(CnGroup),intent(in) :: CgroupStart + integer :: maxSkipsForAtom + #ifdef IS_MPI integer, allocatable, dimension(:) :: c_idents_Row integer, allocatable, dimension(:) :: c_idents_Col @@ -98,12 +96,9 @@ contains nLocal = CnLocal nGlobal = CnGlobal + nGroup = CnGroup thisSim = setThisSim - rcut2 = thisSim%rcut * thisSim%rcut - rcut6 = rcut2 * rcut2 * rcut2 - rlist2 = thisSim%rlist * thisSim%rlist - box = thisSim%box nExcludes_Global = CnGlobalExcludes nExcludes_Local = CnLocalExcludes @@ -166,21 +161,81 @@ contains deallocate(c_idents_Row) endif -#else +#endif + +! We build the local atid's for both mpi and nonmpi do i = 1, nLocal me = getFirstMatchingElement(atypes, "c_ident", c_idents(i)) atid(i) = me enddo -#endif - - do i = 1, nExcludes_Local excludesLocal(1,i) = CexcludesLocal(1,i) excludesLocal(2,i) = CexcludesLocal(2,i) enddo + + maxSkipsForAtom = 0 + do j = 1, nLocal + nSkipsForAtom(j) = 0 +#ifdef IS_MPI + id1 = tagRow(j) +#else + id1 = j +#endif + do i = 1, nExcludes_Local + if (excludesLocal(1,i) .eq. id1 ) then + nSkipsForAtom(j) = nSkipsForAtom(j) + 1 + + if (nSkipsForAtom(j) .gt. maxSkipsForAtom) then + maxSkipsForAtom = nSkipsForAtom(j) + endif + endif + if (excludesLocal(2,i) .eq. id1 ) then + nSkipsForAtom(j) = nSkipsForAtom(j) + 1 + + if (nSkipsForAtom(j) .gt. maxSkipsForAtom) then + maxSkipsForAtom = nSkipsForAtom(j) + endif + endif + end do + enddo + + allocate(skipsForAtom(nLocal, maxSkipsForAtom), stat=alloc_stat) + if (alloc_stat /= 0 ) then + write(*,*) 'Could not allocate skipsForAtom array' + return + endif + + do j = 1, nLocal + nSkipsForAtom(j) = 0 +#ifdef IS_MPI + id1 = tagRow(j) +#else + id1 = j +#endif + do i = 1, nExcludes_Local + if (excludesLocal(1,i) .eq. id1 ) then + nSkipsForAtom(j) = nSkipsForAtom(j) + 1 +#ifdef IS_MPI + id2 = tagColumn(excludesLocal(2,i)) +#else + id2 = excludesLocal(2,i) +#endif + skipsForAtom(j, nSkipsForAtom(j)) = id2 + endif + if (excludesLocal(2, i) .eq. id2 ) then + nSkipsForAtom(j) = nSkipsForAtom(j) + 1 +#ifdef IS_MPI + id2 = tagColumn(excludesLocal(1,i)) +#else + id2 = excludesLocal(1,i) +#endif + skipsForAtom(j, nSkipsForAtom(j)) = id2 + endif + end do + enddo do i = 1, nExcludes_Global excludesGlobal(i) = CexcludesGlobal(i) @@ -188,97 +243,50 @@ contains do i = 1, nGlobal molMemberShipList(i) = CmolMembership(i) - enddo + enddo + + biggestGroupSize = 0 + do i = 1, nGroup + groupStart(i) = CgroupStart(i) + groupSize = 0 + gStart = groupStart(i) + if (i .eq. nGroup) then + gEnd = nLocal + else + gEnd = CgroupStart(i+1) - 1 + endif + do ia = gStart, gEnd + groupList(ia) = CgroupList(ia) + groupSize = groupSize + 1 + enddo + if (groupSize .gt. biggestGroupSize) biggestGroupSize = groupSize + enddo + groupStart(nGroup+1) = nLocal+1 + do i = 1, nLocal + mfact(i) = Cmfact(i) + enddo + if (status == 0) simulation_setup_complete = .true. end subroutine SimulationSetup - subroutine setBox_3d(new_box_size) - real(kind=dp), dimension(3) :: new_box_size + subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic) + real(kind=dp), dimension(3,3) :: cHmat, cHmatInv + integer :: cBoxIsOrthorhombic integer :: smallest, status, i - - thisSim%box = new_box_size - box = thisSim%box - - return - end subroutine setBox_3d - - subroutine setBox_1d(dim, new_box_size) - integer :: dim, status - real(kind=dp) :: new_box_size - thisSim%box(dim) = new_box_size - box(dim) = thisSim%box(dim) - end subroutine setBox_1d - - subroutine setRcut(new_rcut, status) - real(kind = dp) :: new_rcut - integer :: myStatus, status - thisSim%rcut = new_rcut - rcut2 = thisSim%rcut * thisSim%rcut - rcut6 = rcut2 * rcut2 * rcut2 - status = 0 - return - end subroutine setRcut - - function getBox_3d() result(thisBox) - real( kind = dp ), dimension(3) :: thisBox - thisBox = thisSim%box - end function getBox_3d - - function getBox_1d(dim) result(thisBox) - integer, intent(in) :: dim - real( kind = dp ) :: thisBox - thisBox = thisSim%box(dim) - end function getBox_1d - - subroutine getRcut(thisrcut,rc2,rc6,status) - real( kind = dp ), intent(out) :: thisrcut - real( kind = dp ), intent(out), optional :: rc2 - real( kind = dp ), intent(out), optional :: rc6 - integer, optional :: status - - if (present(status)) status = 0 + Hmat = cHmat + HmatInv = cHmatInv + if (cBoxIsOrthorhombic .eq. 0 ) then + boxIsOrthorhombic = .false. + else + boxIsOrthorhombic = .true. + endif - if (.not.simulation_setup_complete ) then - if (present(status)) status = -1 - return - end if - - thisrcut = thisSim%rcut - if(present(rc2)) rc2 = rcut2 - if(present(rc6)) rc6 = rcut6 - end subroutine getRcut - - subroutine getRlist(thisrlist,rl2,status) - real( kind = dp ), intent(out) :: thisrlist - real( kind = dp ), intent(out), optional :: rl2 + return + end subroutine setBox - integer, optional :: status - - if (present(status)) status = 0 - - if (.not.simulation_setup_complete ) then - if (present(status)) status = -1 - return - end if - - thisrlist = thisSim%rlist - if(present(rl2)) rl2 = rlist2 - end subroutine getRlist - - function getRrf() result(rrf) - real( kind = dp ) :: rrf - rrf = thisSim%rrf - write(*,*) 'getRrf = ', rrf, thisSim%rrf - end function getRrf - - function getRt() result(rt) - real( kind = dp ) :: rt - rt = thisSim%rt - end function getRt - function getDielect() result(dielect) real( kind = dp ) :: dielect dielect = thisSim%dielect @@ -299,6 +307,11 @@ contains doesit = thisSim%SIM_uses_sticky end function SimUsesSticky + function SimUsesCharges() result(doesit) + logical :: doesit + doesit = thisSim%SIM_uses_charges + end function SimUsesCharges + function SimUsesDipoles() result(doesit) logical :: doesit doesit = thisSim%SIM_uses_dipoles @@ -342,6 +355,12 @@ contains thisStat = 0 call FreeSimGlobals() + + allocate(nSkipsForAtom(nLocal), stat=alloc_stat) + if (alloc_stat /= 0 ) then + thisStat = -1 + return + endif allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat) if (alloc_stat /= 0 ) then @@ -360,6 +379,24 @@ contains thisStat = -1 return endif + + allocate(groupStart(nGroup+1), stat=alloc_stat) + if (alloc_stat /= 0 ) then + thisStat = -1 + return + endif + + allocate(groupList(nLocal), stat=alloc_stat) + if (alloc_stat /= 0 ) then + thisStat = -1 + return + endif + + allocate(mfact(nLocal), stat=alloc_stat) + if (alloc_stat /= 0 ) then + thisStat = -1 + return + endif end subroutine InitializeSimGlobals @@ -367,16 +404,20 @@ contains !We free in the opposite order in which we allocate in. + if (allocated(skipsForAtom)) deallocate(skipsForAtom) + if (allocated(mfact)) deallocate(mfact) + if (allocated(groupList)) deallocate(groupList) + if (allocated(groupStart)) deallocate(groupStart) if (allocated(molMembershipList)) deallocate(molMembershipList) if (allocated(excludesGlobal)) deallocate(excludesGlobal) if (allocated(excludesLocal)) deallocate(excludesLocal) - + end subroutine FreeSimGlobals - + pure function getNlocal() result(n) integer :: n n = nLocal end function getNlocal - + end module simulation