Scrambler  1
hyperbolic_control.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 !    hyperbolic_control.f90 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 
00031 
00035 
00036 
00039 
00040 !============================================================================================
00041 ! Module Name:        SolverControl
00042 ! Module File:        solver_control.f90
00043 ! Purpose:              Provide an interface for the rest of the code to use when 
00044 !                       initializing and using numerical solvers.  Handles the 
00045 !                       selection of a scheme and the reading in of appropriate data.
00046 ! Public Methods:    SolverInit(), advance()
00047 ! Created:            20100705 by Brandon D. Shroyer
00048 !============================================================================================
00049 MODULE HyperbolicControl
00050 
00051   USE DataDeclarations
00052   USE HyperbolicDeclarations
00053   USE PhysicsDeclarations
00054   USE SweepScheme
00055   USE MUSCLScheme
00056   USE Timing
00057   IMPLICIT NONE
00058   PRIVATE
00059 
00060   PUBLIC HyperbolicInit, advance, PRofileAdvance, AdvanceStackSize
00061 
00062 CONTAINS
00063 
00069 
00070 
00076   SUBROUTINE HyperbolicInit()
00077     INTEGER :: iErr
00078     ALLOCATE(maxspeed(0:MaxLevel), maxwavespeed(0:Maxlevel), maxsolverspeed(0:MaxLevel), AdvanceStencil(0:MaxLevel), t_startadvance(0:MaxLevel))
00079 
00080 
00081     ALLOCATE(WorkDoneByLevel(-1:MaxLevel))
00082     ALLOCATE(WorkDoneByGrid(-1:MaxLevel))
00083     ALLOCATE(InternalCellUpdates(0:MaxLevel))
00084     ALLOCATE(CellUpdates(0:MaxLevel))
00085 
00086     WorkDoneByLevel=0d0
00087     WorkDoneByGrid=0d0
00088     InternalCellUpdates=0
00089     CellUpdates=0
00090 
00091     ! Open data file (quit if not successful).
00092     OPEN(UNIT=SOLVER_DATA_HANDLE,FILE=SOLVER_DATA_FILE,IOSTAT=iErr)
00093 
00094     IF (iErr /= 0) THEN
00095        PRINT *, "SolverInit() error: unable to open ", SOLVER_DATA_FILE, "."
00096        STOP
00097     END IF
00098 
00099     ! Read in scheme value.
00100     READ(SOLVER_DATA_HANDLE,NML=SolverData,IOStat=iErr)
00101 
00102     IF(iErr/=0) THEN
00103        PRINT *, "SolverInit() error:  unable to read ", SOLVER_DATA_FILE, "."
00104        STOP
00105     END IF
00106 
00107 
00108     IF ((iSolver == 1 .OR. iSolver == 3) .AND. (lMHD .OR. (iEOS == EOS_ISOTHERMAL))) THEN
00109        PRINT*,'method(1) = 1 and iSolver = ', iSolver, 'is only supported for adiabatic hydro.  ', &
00110               'Please verify that iCooling /= 4 and that lMHD = False '
00111     END IF
00112     IF (iSolver== 0) THEN
00113        iSolver=2
00114        IF (lMHD) iSolver = iSolver + 4
00115        IF (iEOS == EOS_ISOTHERMAL) iSolver = iSolver + 2
00116     END IF
00117 
00118     ! Select the appropriate scheme option--currently there are only two, so
00119     ! there's no need for a select-case statement.
00120     IF (iScheme == SWEEP_SCHEME_ID) THEN
00121        CALL SweepReadDomainData()
00122 
00123     ELSE IF (iScheme == MUSCL_SCHEME_ID) THEN
00124        CALL MUSCLInit()
00125 
00126     ELSE
00127        PRINT "('SolverInit() error: invalid scheme option ', i3, '.')", iScheme
00128        STOP
00129 
00130     END IF
00131 
00132     CLOSE(SOLVER_DATA_HANDLE, IOStat=iErr)
00133 
00134     IF (iErr /= 0) THEN
00135        PRINT *, "SolverInit() error:  unable to close ", SOLVER_DATA_FILE, "."
00136     END IF
00137 
00138 
00139     ALLOCATE(NodeCompleted(0:MaxLevel))
00140     ALLOCATE(tused_this_grid(0:MaxLevel))
00141   END SUBROUTINE HyperbolicInit
00142 
00143 
00144   FUNCTION AdvanceStackSize(n)
00145      INTEGER :: n
00146      INTEGER :: AdvanceStackSize
00147      SELECT CASE(iScheme)
00148      CASE(SWEEP_SCHEME_ID)
00149         AdvanceStackSize=SweepAdvanceStackSize(n)
00150      CASE DEFAULT
00151         PRINT*, 'error - ischeme not recognized'
00152      END SELECT
00153 
00154   END FUNCTION AdvanceStackSize
00155 
00161   SUBROUTINE advance(Info, lCompleteOpt)
00162     ! Interface declarations
00163     TYPE (InfoDef) :: Info
00164     LOGICAL :: partialOK, lComplete
00165     INTEGER :: bc, mB(3,2), i
00166     LOGICAL, OPTIONAL :: lCompleteOpt
00167     IF (Present(lCompleteOpt)) THEN
00168        lComplete=lCompleteOpt
00169     ELSE
00170        lComplete=.false.
00171     END IF
00172     ! choose which scheme
00173     !    bc(1:nDim) = levels(level)%gmbc(levels(level)%step)-hyperbolic_mbc
00174     !    dom_range= reshape((/1, 1, 1, Info%mX(1),Info%mX(2), Info%mX(3)/),(/3,2/))
00175     !    dom_range=expand(1,bc(1), expand(2,bc(2), expand(3,bc(3), dom_range)))
00176     ! Number of ghost zones that need to be updated by hyperbolic advance
00177 !    if (levels(Info%level)%dt == 0) THEN
00178 !       write(*,*) "shouldn't be here"
00179 !       write(*,*) bc
00180 !       STOP
00181 !    end if
00182     
00183     bc=levels(Info%level)%ambc(levels(Info%level)%step) 
00184     mB=1
00185     mB(1:nDim,2)=Info%mX(1:nDim)+bc
00186     mB(1:nDim,1)=1-bc
00187 
00188 !    DO i=1, 100
00189 !       Info%q(mb(1,1):mb(1,2), mb(2,1):mb(2,2), mb(3,1):mb(3,2),1:NrHydroVars)=2d0*&
00190 !            Info%q(mb(1,1):mb(1,2), mb(2,1):mb(2,2), mb(3,1):mb(3,2),1:NrHydroVars)
00191 !       Info%q(mb(1,1):mb(1,2), mb(2,1):mb(2,2), mb(3,1):mb(3,2),1:NrHydroVars)=.5d0*&
00192 !            Info%q(mb(1,1):mb(1,2), mb(2,1):mb(2,2), mb(3,1):mb(3,2),1:NrHydroVars)
00193 !       Info%q(mb(1,1):mb(1,2), mb(2,1):mb(2,2), mb(3,1):mb(3,2),1:NrHydroVars)=1.001*&
00194 !            Info%q(mb(1,1):mb(1,2), mb(2,1):mb(2,2), mb(3,1):mb(3,2),1:NrHydroVars)
00195 !    END DO
00196 !    RETURN
00197 !    IF (MPI_ID == 0) write(*,*) 'advancing', mB
00198     SELECT CASE(iScheme)
00199     CASE(MUSCL_SCHEME_ID)
00200        CALL MUSCLAdvance(Info,mB, lComplete)
00201     CASE(SWEEP_SCHEME_ID)
00202        CALL sweepAdvance(Info,mB, lComplete)
00203     CASE DEFAULT
00204        PRINT "('advance() error: invalid scheme option ', i3, '.')", iScheme
00205        STOP
00206 
00207     END SELECT
00208 !    Info%OldCostPerCell=Info%CostPerCell
00209 !    Info%CostPerCell= 1
00210   END SUBROUTINE advance
00211 
00212 
00213   SUBROUTINE ProfileAdvance
00214      USE ModuleControl
00215      USE CommonFunctions
00216      !    USE IOBOV
00217      INTEGER :: i,j,k,l,s, ierr, m
00218      INTEGER :: ncase,nvar,nmax,nk,nsize, npoints
00219      INTEGER, DIMENSION(3,2) :: mGlobal
00220      REAL(KIND=qPREC) :: tstart, temp, temp2
00221      TYPE(InfoDef), POINTER :: Info
00222      REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: Times,sls,errls
00223      INTEGER, DIMENSION(:), ALLOCATABLE :: points
00224      REAL(KIND=qPREC), DIMENSION(:,:), ALLOCATABLE :: LeastSquares
00225      REAL(KIND=qPREC) :: AdvanceCoeffsTemp(8), err
00226      CHARACTER(LEN=40) :: FileName
00227      INTEGER :: filehandle=1324, tempithreaded, curr_progress, old_progress
00228      REAL(KIND=qPrec), DIMENSION(:), ALLOCATABLE :: MeanAdvanceCoeffs
00229      LOGICAL :: lFound
00230      NAMELIST /ProfileData/ AdvanceCoeffsTemp     ! This namelist will probably grow over time.   
00231      levels(0)%step=levels(0)%steps
00232      !     levels(0)%gmbc(1)=hyperbolic_mbc
00233      !     levels(0)%dx=1
00234      levels(0)%dt=1e-10
00235      tempiThreaded=iThreaded
00236      iThreaded = NON_THREADED
00237      ! This section computes the least square from sampled mpi_wtime
00238      !     ncase=ndim**(6-ndim);nvar=2**nDim;nmax=nDim**(9-2*nDim);nk=27**(nDim-2);nsize=6-ndim
00239 
00240 
00241      IF (nDim == 1) THEN
00242         npoints=12
00243         nvar=2
00244         ALLOCATE(AdvanceGridTimes(npoints,1,1))
00245      ELSE IF (nDim == 2) THEN
00246         npoints=11
00247         ALLOCATE(AdvanceGridTimes(npoints,npoints,1))
00248         nvar=2
00249      ELSE
00250         npoints=6
00251         ALLOCATE(AdvanceGridTimes(npoints,npoints,npoints))
00252         nvar=1
00253      END IF
00254      ncase=min(4,npoints)**nDim
00255 
00256      ALLOCATE(AdvanceCoeffs(nvar), MeanAdvanceCoeffs(nvar)) !Advance Coefficients
00257      OPEN(UNIT=PROFILE_DATA_HANDLE, IOSTAT=ierr, file=PROFILE_DATA_FILE, status="old")
00258      lFound=.false.
00259      IF (ierr==0) THEN !Found the file
00260         IF (MPI_ID == 0) PRINT *,'Found existing profile.data file.'
00261         lFound=.true.
00262         READ(PROFILE_DATA_HANDLE, NML=ProfileData)
00263         AdvanceCoeffs(1)=AdvanceCoeffsTemp(1)
00264         AdvanceCoeffs(2:size(AdvanceCoeffs))=AdvanceCoeffsTemp(2:size(AdvanceCoeffs))*AdvanceCoeffs(1)
00265         READ(PROFILE_DATA_HANDLE, '(10E25.16)') AdvanceGridTimes
00266         CLOSE(PROFILE_DATA_HANDLE)
00267         IF (.NOT. lSkipProfile) THEN
00268            ! Test profile.data quickly for a small 32**nDim grid
00269            mGlobal=1
00270            ! Get Advance Times (solution) for sample points
00271            mGlobal(1:nDim,2)=2**(5)
00272            CALL InitInfo(Info, 0, mGlobal)
00273            IF (MaintainAuxArrays) Info%aux(:,:,:,:)=0
00274            Info%q(:,:,:,:)=0d0
00275            Info%q(:,:,:,1)=1d0
00276            IF (iE .ne. 0d0) Info%q(:,:,:,iE)=1d0
00277            CALL sweepAdvance(Info,mGlobal,lComplete=.true., lProfile_opt=.true.)
00278            tstart=MPI_Wtime()
00279            CALL sweepAdvance(Info,mGlobal,lComplete=.true., lProfile_opt=.true.)
00280            tstart=MPI_Wtime()-tstart
00281            CALL DestroyInfo(Info)
00282            err=ABS(tstart-AdvanceCost(mGlobal(:,2)))/tstart
00283            CALL MPI_ALLREDUCE(MPI_IN_PLACE, err, 1, MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD, iErr)
00284            IF (err > .1d0) THEN
00285               lFound=.false.
00286            END IF
00287            IF (MPI_ID == 0 .AND. .NOT. lFound) THEN
00288               write(*,*) 'profile.data is off by ', nint(err*100d0), '%. Generating a new profile'
00289               write(*,*) 'If you wish to skip re-profiling - consider adding lSkipProfile=.true. to your global.data file'
00290            END IF
00291         END IF
00292      END IF
00293 
00294      !       IF (MPI_ID == 0) CALL WriteBOV3DScalar('advancecoeffs', (/0d0,0d0,0d0/)-half,(/npoints, npoints, npoints/)-half, 0d0, AdvanceGridTimes, 'advance_times')
00295      IF (.NOT. lfound) THEN
00296         ALLOCATE(Times(1:ncase)) ! Array to store solution for fit points
00297         ALLOCATE(LeastSquares(ncase,nvar))
00298         AdvanceCoeffs=0
00299         mGlobal=1
00300         ! Get Advance Times (solution) for sample points
00301         mGlobal(1:nDim,2)=2**(5)
00302         IF (.NOT. lSkipProfile) THEN
00303            CALL InitInfo(Info, 0, mGlobal)
00304            IF (MaintainAuxArrays) Info%aux(:,:,:,:)=0
00305            Info%q(:,:,:,:)=0d0
00306            Info%q(:,:,:,1)=1d0
00307            IF (iE .ne. 0d0) Info%q(:,:,:,iE)=1d0
00308            CALL sweepAdvance(Info,mGlobal,lComplete=.true.,lProfile_opt=.true.)
00309            CALL DestroyInfo(Info)
00310         END IF
00311         CALL MPI_BARRIER(MPI_COMM_WORLD, ierr) !want these to happen at the same time
00312 
00313         IF (lSkipProfile) THEN
00314            DO i=1,npoints
00315               DO j=1, npoints        
00316                  DO k=1, npoints
00317                     IF (nDim == 1 .AND. j > 1) CYCLE
00318                     IF (nDim <= 2 .AND. k > 1) CYCLE
00319                     AdvanceGridTimes(i,j,k)=2**(i+j+k-3)
00320                  END DO
00321               END DO
00322            END DO
00323         ELSE
00324            IF (MPI_ID == 0)  write(*,*) 'profiling...'
00325            s=1
00326            l=0
00327            curr_progress=0
00328            DO i=1, npoints
00329               DO j=1, npoints        
00330                  DO k=1, npoints
00331                     IF (nDim == 1 .AND. j > 1) CYCLE
00332                     IF (nDim <= 2 .AND. k > 1) CYCLE
00333                     mglobal(1:3,2)=2**((/i,j,k/)-1)
00334                     IF (nDim <= 2) mGlobal(3,2)=1
00335                     IF (nDim <= 1) mGlobal(2,2)=1
00336                     CALL InitInfo(Info, 0, mGlobal)
00337                     Info%q(:,:,:,:)=0d0
00338                     Info%q(:,:,:,1)=1d0
00339                     IF (iE .ne. 0d0) Info%q(:,:,:,iE)=1d0
00340                     IF (MaintainAuxArrays) Info%aux=0d0
00341                     tstart=MPI_Wtime()
00342                     CALL sweepAdvance(Info,mGlobal,lComplete=.true.,lProfile_opt=.true.)
00343                     AdvanceGridTimes(i,j,k)=MPI_Wtime()-tstart
00344                     !                 IF (i > npoints-4 .AND. (nDim < 2 .OR. j > npoints-4) .AND. (nDim < 3 .OR. k > npoints-4)) THEN
00345                     !                    Times(s)=AdvanceGridTimes(i,j,k)
00346                     !                    LeastSquares(s,:)=(/product(Info%mX), 1/)
00347                     !                    s=s+1
00348                     !                 END IF
00349                     CALL DestroyInfo(Info)
00350                     l=l+1
00351                     curr_progress=nint(40d0*real(l)/real(npoints**nDim))
00352                     IF (MPI_ID == 0) write(*,'(1A1,A1,I3,A2,40A1,A1,$)')  char(13), ' ', nint(100d0*real(curr_progress)/40d0), '%|', ('|', m=1, curr_progress), (' ', m=curr_progress+1,40), '|'
00353                  END DO
00354               END DO
00355            END DO
00356            IF (MPI_ID == 0) write(*,*)
00357         END IF
00358 
00359         !Call the least square fitting function, 
00360         !sls stores the solution vector in an order of (xy,x,y,1) or (xyz,xy,yz,zx,x,y,z,1)
00361         !errls stores the relative error for the ncase sampled data
00362 !        CALL LSSolver(LeastSquares,Times,AdvanceCoeffs) ! 
00363 
00364         AdvanceCoeffs(1)=AdvanceGridTimes(npoints,merge(npoints,1,nDim >= 2) ,merge(npoints,1,nDim == 3))/(2d0**(npoints-1))**ndim
00365 
00366         DEALLOCATE(LeastSquares,Times) ! Array to store solution for fit points
00367 
00368         mGlobal=1
00369 
00370         IF (MPI_NP > 1) THEN
00371            CALL MPI_ALLREDUCE(AdvanceCoeffs, MeanAdvanceCoeffs, size(AdvanceCoeffs), MPI_DOUBLE_PRECISION, MPI_SUM, MPI_COMM_WORLD, iErr)
00372            CALL MPI_ALLREDUCE(MPI_IN_PLACE, AdvanceGridTimes, size(AdvanceGridTimes), MPI_DOUBLE_PRECISION, MPI_SUM, MPI_COMM_WORLD, iErr)
00373            AdvanceGridTimes=AdvanceGridTimes/REAL(MPI_NP)
00374            MeanAdvanceCoeffs=MeanAdvanceCoeffs/REAL(MPI_NP)
00375            AdvanceCoeffs=MeanAdvanceCoeffs
00376         END IF
00377         AdvanceGridTimes=log(AdvanceGridTimes)
00378 
00379         IF (MPI_ID == 0 .AND. .NOT. lSkipProfile) THEN 
00380            AdvanceCoeffsTemp=0d0
00381            AdvanceCoeffsTemp(1)=AdvanceCoeffs(1)
00382            AdvanceCoeffsTemp(2:size(AdvanceCoeffs))=AdvanceCoeffs(2:size(AdvanceCoeffs))/AdvanceCoeffs(1)
00383            open(UNIT=PROFILE_DATA_HANDLE, file=PROFILE_DATA_FILE, status="unknown")
00384            WRITE(PROFILE_DATA_HANDLE, NML=ProfileData)
00385            WRITE(PROFILE_DATA_HANDLE, '(10E25.16)') AdvanceGridTimes
00386            CLOSE(PROFILE_DATA_HANDLE)
00387         END IF
00388 
00389      END IF
00390      IF (1 == 2) THEN
00391         CALL MPI_BARRIER(MPI_COMM_WORLD, iErr)
00392         write(*,*) "AdvanceCoeffs = ", AdvanceCoeffs
00393         DO i=0,10,1
00394            DO j=0,10,1
00395               DO k=1,1,1
00396                  IF (nDim == 2 .AND. k /= 1) CYCLE
00397                  mGlobal(1:3,2)=(/2**i,2**j,k/)
00398                  !                mGlobal(1:3,2)=(/i,j,k/)
00399                  CALL InitInfo(Info, 0, mGlobal)
00400                  CALL GridInit(Info)
00401                  tstart=MPI_Wtime()
00402                  CALL sweepAdvance(Info,mGlobal,lComplete=.true.,lProfile_opt=.true.)
00403                  temp=MPI_Wtime()-tstart
00404                  temp2=SimpleAdvanceCost(mGlobal(:,2))
00405                  temp2=AdvanceCost(mGlobal(:,2))
00406                  !                write(*,'(2I6,2E13.5,F13.5)') 2**i,2**j,temp,temp2, 2d0*((temp-temp2)/(temp+temp2))*100
00407                  write(*,'(2I6,2E13.5,F13.5)') 2**i,2**j, 2d0*((temp-temp2)/(temp+temp2))*100
00408                  !                write(*,*) i,j,temp2
00409                  !                write(*,*) i,j,temp2 !temp2/(Product(REAL(mglobal(1:nDim,2)))*AdvanceCoeffs(1))
00410                  CALL DestroyInfo(Info)
00411               END DO
00412            END DO
00413            write(*,*) 
00414         END DO
00415      END IF
00416      iThreaded = TempiThreaded
00417 
00418 
00419      AdvanceTimer%Accumulator(:)=0
00420      AdvancePredictor%Accumulator(:)=0
00421      MySpeedFactor=1d0
00422      !    CALL MPI_BARRIER(MPI_COMM_WORLD, ierr)
00423      !    STOP
00424   END SUBROUTINE ProfileAdvance
00425 
00426 
00427 END MODULE HyperbolicControl
00428 
00429 
 All Classes Files Functions Variables