Scrambler  1
FillingFraction/problem.f90
Go to the documentation of this file.
00001 !#########################################################################
00002 !               
00003 !    Copyright (C) 2003-2012 Department of Physics and Astronomy,
00004 !                            University of Rochester,
00005 !                            Rochester, NY
00006 !
00007 !    problem.f90 of module FillingFraction is part of AstroBEAR.
00008 !
00009 !    AstroBEAR is free software: you can redistribute it and/or modify    
00010 !    it under the terms of the GNU General Public License as published by 
00011 !    the Free Software Foundation, either version 3 of the License, or    
00012 !    (at your option) any later version.
00013 !
00014 !    AstroBEAR is distributed in the hope that it will be useful, 
00015 !    but WITHOUT ANY WARRANTY; without even the implied warranty of
00016 !    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00017 !    GNU General Public License for more details.
00018 !
00019 !    You should have received a copy of the GNU General Public License
00020 !    along with AstroBEAR.  If not, see <http://www.gnu.org/licenses/>.
00021 !
00022 !#########################################################################
00025 
00028 
00032 
00035 MODULE Problem
00036    USE DataDeclarations
00037    USE SplitRegions
00038    USE Shapes
00039    USE EOS
00040    USE RiemannSolvers
00041    USE Ambients
00042    USE TreeDeclarations
00043    IMPLICIT NONE
00044    SAVE
00045 
00046    PUBLIC ProblemModuleInit, ProblemGridInit, &
00047         ProblemBeforeStep, ProblemAfterStep, ProblemSetErrFlag, ProblemBeforeGlobalStep
00048    PRIVATE
00049    REAL(KIND=qPREC), DIMENSION(0:MaxDepth) :: filling_fractions
00050    REAL(KIND=qPREC) :: vel(3), L
00051    REAL(8) :: ProbtStart
00052    REAL(8) :: ProbTFinal
00053    INTEGER(8) :: nRootSteps, iRootSteps, StartStep, CellUpdatesByLevel(0:MaxDepth), TotalCellUpdatesByLevel(0:MaxDepth)
00054    INTEGER :: branching_ratio(0:MaxDepth)=1
00055    REAL(KIND=qPREC) :: CellCostByLevel(0:MaxDepth), TotalCellCostByLevel(0:MaxDepth), dy, dz
00056    REAL(KIND=qPREC), ALLOCATABLE, DIMENSION(:,:,:,:) :: GxByLevel
00057    !  REAL(KIND=qPREC), ALLOCATABLE, DIMENSION(:,:) :: WorkLoadByStepByLevel
00058    INTEGER, ALLOCATABLE, DIMENSION(:,:) :: CellsByStepByLevel
00059    REAL(KIND=qPREC), ALLOCATABLE, DIMENSION(:) :: WorkLoadByRootStep
00060    INTEGER, ALLOCATABLE, DIMENSION(:) :: CellsByRootStep, indices
00061    LOGICAL :: lPlaceCentered=.true.
00062 CONTAINS
00063 
00065    SUBROUTINE ProblemModuleInit()      
00066       INTEGER :: i
00067       REAL(KIND=qPREC) :: rs_old, rs_new
00068       NAMELIST /ProblemData/ filling_fractions, vel, branching_ratio, lPlaceCentered
00069       OPEN(UNIT=PROBLEM_DATA_HANDLE, FILE='problem.data', STATUS="OLD")
00070       READ(PROBLEM_DATA_HANDLE,NML=ProblemData)
00071       CLOSE(PROBLEM_DATA_HANDLE)
00072       filling_fractions=filling_fractions(0)
00073 
00074       L=.5*(GxBounds(1,2)-GxBounds(1,1))
00075       levels(0)%qTolerance=filling_fractions(0)**(-1d0/REAL(nDim))/L*levels(0)%dx
00076       DO i=1, MaxLevel
00077          levels(i)%qTolerance=.5d0*levels(i-1)%qTolerance*filling_fractions(i)**(-1d0/REAL(nDim))
00078       END DO
00079       !      write(*,*) levels(:)%qTolerance
00080       iRootSteps=0
00081       dy=(GxBounds(2,2)-GxBounds(2,1))/(GxBounds(2,2)-GxBounds(2,1))
00082       IF (nDim == 2) THEN
00083          dz=0
00084       ELSE 
00085          dz=(GxBounds(1,2)-GxBounds(1,1))/(GxBounds(3,2)-GxBounds(3,1))
00086       END IF
00087 
00088       ALLOCATE(GxByLevel(0:MaxLevel,product(branching_ratio(0:MaxLevel-1)),3,2))
00089       GxByLevel=1
00090       ALLOCATE(indices(0:MaxLevel))
00091       indices=0
00092      
00093       CALL CreateClumpLets(half*SUM(GxBounds,2), 0, half*(GxBounds(1,2)-GxBounds(1,1)))
00094    END SUBROUTINE ProblemModuleInit
00095 
00096    RECURSIVE SUBROUTINE CreateClumplets(pos,n,hw)
00097       REAL(KIND=qPREC), DIMENSION(3) :: pos
00098       INTEGER :: n,i
00099       REAL(KIND=qPREC), DIMENSION(:,:), ALLOCATABLE :: childpos
00100       REAL(KIND=qPREC) :: hw, hwc
00101       indices(n)=indices(n)+1
00102       GxByLevel(n,indices(n),1:nDim,1)=pos(1:nDim)-hw
00103       GxByLevel(n,indices(n),1:nDim,2)=pos(1:nDim)+hw
00104       IF (MPI_ID == 0) write(*,'(6F12.5)') GxByLevel(n,indices(n),1:nDim,1:2), hw
00105       IF (n < MaxLevel) THEN
00106          ALLOCATE(childpos(branching_ratio(n),3))
00107          hwc=hw*(filling_fractions(n)/real(branching_ratio(n)))**(1d0/REAL(nDim))
00108          IF (branching_ratio(n) == 1 .AND. lPlaceCentered) THEN
00109             childpos(1,:)=pos
00110          ELSE
00111             CALL GetChildPos(branching_ratio(n),childpos,hwc/hw*1.1) !REAL(nDim)))
00112 !         write(*,*) childpos, hwc
00113          END IF
00114          DO i=1,branching_ratio(n)
00115             CALL CreateClumpLets(pos+hw*childpos(i,:), n+1, hwc)
00116          END DO
00117          DEALLOCATE(childpos)
00118       END IF
00119    END SUBROUTINE CreateClumplets
00120 
00121 
00122    SUBROUTINE GetChildPos(n_children,childpos,r)
00123       INTEGER :: n_children
00124       REAL(KIND=qPREC), DIMENSION(:,:) :: childpos
00125       REAL(KIND=qPREC) :: r
00126       INTEGER, PARAMETER :: nIterations=100000
00127       LOGICAL :: lSuccess
00128       INTEGER :: j,k,l
00129       REAL :: a
00130       !create childpositions where each child is at least 2r apart and within (1-r) of the center...
00131       childpos=0
00132       DO j=1,nIterations
00133          lSuccess=.true.
00134          DO k=1,n_children
00135             CALL Random(a)
00136             childpos(k,1)=a
00137             CALL Random(a)
00138             childpos(k,2)=a
00139             IF (nDim == 3) THEN
00140                CALL Random(a)
00141                childpos(k,3)=a
00142             END IF
00143             Childpos(k,1:nDim)=2d0*(childpos(k,1:nDim)-.5d0)*(1d0-r) !rescale to be within a random cube with half width (1-r)                      
00144             DO l=1,k-1
00145                IF (ALL(childpos(k,1:nDim)-childpos(l,1:nDim) < 2d0*r)) THEN
00146                   lSuccess=.false.
00147                   EXIT
00148                END IF
00149             END DO
00150             IF (.Not. lSuccess) EXIT
00151          END DO
00152          IF (lSuccess) EXIT
00153       END DO
00154       IF (.NOT. lSuccess) THEN
00155          IF (MPI_ID == 0) write(*,*) 'failed to place', n_children, 'clumplets of radius ',r
00156 
00157          STOP
00158       END IF
00159 
00160    END SUBROUTINE GetChildPos
00161 
00164    SUBROUTINE ProblemGridInit(Info)
00165       TYPE(InfoDef) :: Info
00166       INTEGER :: i,j,k, rmbc, zrmbc
00167       REAL(KIND=qPREC) :: x, y, z, dx, rho, r
00168       Info%q=0
00169       Info%q(:,:,:,1)=Info%level+1
00170       DO i=1,nDim
00171          Info%q(:,:,:,imom(i))=(Info%level+1)*vel(i)
00172       END DO
00173       Info%q(:,:,:,iE)=1d0+(Info%level+1)*sum(vel**2)/2d0
00174    END SUBROUTINE ProblemGridInit
00175 
00176 
00177    SUBROUTINE ProblemBeforeGlobalStep(n)
00178       USE Timing
00179       REAL(8) :: temp
00180       INTEGER :: n, i, LocalNodes, ExtNodes, TotalNodes, MyCells, TotalCells, TotalExtNodes, iErr
00181       REAL(8) :: AdvanceTime(9)
00182       INTEGER, ALLOCATABLE, DIMENSION(:,:) :: TempMeanCellsByStepByLevel, MaxCellsByStepByLevel, CellsByRootStep
00183       REAL(KIND=qPREC), ALLOCATABLE, DIMENSION(:,:) :: MaxWorkLoadByStepByLevel, MeanWorkLoadByStepByLevel, MeanCellsByStepByLevel, sigmacellsbystepbylevel, sigmaworkbystepbylevel, workbyrootstep
00184       INTEGER :: nFinestSteps
00185       REAL(KIND=qPREC) :: mystats(3), totalstats(3)
00186 
00187       !      write(*,*) n, iRootSteps, StartStep, nRootSteps        
00188 !      IF (n == 0) THEN
00189          !         write(*,*) levels(n)%CurrentLevelStep
00190 !         IF (levels(0)%CurrentLevelStep==StartStep) THEN
00191 !            ProbtStart= MPI_WTime()
00192 !            CellUpdatesByLevel=0
00193 !            AdvanceTimer%Accumulator(:)=0d0
00194 !            Timers(iBarrier)%Accumulator(:)=0d0
00195 !            Timers(iAmr)%Accumulator(:)=0d0
00196 !         END IF
00197 !      END IF
00198    END SUBROUTINE ProblemBeforeGlobalStep
00199 
00202    SUBROUTINE ProblemBeforeStep(Info)
00203       TYPE(InfoDef) :: Info
00204       INTEGER :: mx(3)
00205       REAL(KIND=qPREC) :: advcost
00206  !     IF (levels(0)%CurrentLevelStep < StartStep) RETURN
00207 !!      mx=Info%mX
00208 !!      CellUpdatesByLevel(Info%level)=CellUpdatesByLevel(Info%level)+PRODUCT(Info%mX)
00209 
00210    END SUBROUTINE ProblemBeforeStep
00211 
00214    SUBROUTINE ProblemAfterStep(Info)
00215       TYPE(InfoDef) :: Info
00216    END SUBROUTINE ProblemAfterStep
00217 
00220    SUBROUTINE ProblemSetErrFlag(Info)
00221       TYPE(InfoDef) :: Info
00222       INTEGER :: sample_res(3), nOverlaps,n,i
00223       INTEGER, DIMENSION(3,2) :: mS    
00224       INTEGER, POINTER, DIMENSION(:,:,:) :: mSs
00225       REAL(KIND=qPREC), DIMENSION(3) :: offset
00226       REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: offsets
00227       DO i=1,indices(Info%level+1)
00228          CALL CalcPhysicalOverlaps(Info, GxByLevel(Info%level+1,i,:,:)+spread(vel,2,2)*levels(Info%level)%tnow, mSs, nOverlaps, offsets, iEVERYWHERE, lHydroPeriodic,0)            
00229          IF (nOverlaps > 0) THEN
00230             DO n=1,nOverlaps
00231                mS=mSs(n,:,:)
00232                Info%ErrFlag(mS(1,1):mS(1,2), mS(2,1):mS(2,2), mS(3,1):mS(3,2))=1
00233             END DO
00234             DEALLOCATE(mSs, offsets)
00235          END IF
00236       END DO
00237    END SUBROUTINE ProblemSetErrFlag
00238 
00239 END MODULE Problem
00240 
 All Classes Files Functions Variables