Scrambler  1
MUSCL_scheme.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 !    MUSCL_scheme.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 
00032 
00033 !===============================================================================
00034 ! Module Name:          MUSCLScheme
00035 ! Module File:          MUSCL_scheme.f90
00036 ! Purpose:              Implements a numerical scheme for solving 
00037 !                       hydrodynamic equations.
00038 ! Public Methods:       ReadMUSCLDomainData(), MUSCLadvance()
00039 ! Created:              by Jonathan Carroll-Nellenback and Baowei Liu 
00040 !===============================================================================
00041 
00046 MODULE MUSCLScheme
00047    USE GlobalDeclarations
00048    USE DataDeclarations
00049    USE HyperbolicDeclarations
00050    USE SourceControl
00051    USE SourceDeclarations
00052    USE PhysicsDeclarations
00053    USE EOS
00054    USE ModuleControl
00055    USE RiemannSolvers
00056    USE DataInfoOps
00057    USE TreeDeclarations
00058    IMPLICIT NONE
00059    PRIVATE
00060    SAVE 
00061    PUBLIC MUSCLInit, MUSCLadvance
00062    INTEGER :: nFields,nSteps
00063    INTEGER, DIMENSION(:,:), ALLOCATABLE :: wdi, fdi, dim_orders
00064    LOGICAL :: lUseLimiter=.true.
00065 CONTAINS
00066    
00067    SUBROUTINE MUSCLInit()
00068       NAMELIST/MUSCLData/ lUseLimiter
00069       INTEGER :: i,d 
00070       READ(SOLVER_DATA_HANDLE,NML=MUSCLData)
00071       hyperbolic_mbc=2
00072 
00073       IF (nDim == 1) THEN      
00074          nSteps=1
00075          ALLOCATE(dim_orders(nDim,nSteps))
00076          dim_orders=reshape((/1/),(/nDim,nSteps/))
00077          nFields=3
00078          ALLOCATE(wdi(nDim,nFields+NrTracerVars), fdi(nDim,nFields+NrTracerVars))
00079          wdi(1,1:nFields)=(/1,3,2/)
00080          fdi(1,1:nFields)=(/1,3,2/)
00081       ELSEIF (nDim == 2) THEN
00082          nSteps=2
00083          ALLOCATE(dim_orders(nDim,nSteps))
00084          dim_orders=reshape((/1,2, 2,1/),(/nDim,nSteps/))
00085 
00086          nFields=4
00087          ALLOCATE(wdi(nDim,nFields+NrTracerVars), fdi(nDim,nFields+NrTracerVars))
00088          wdi(1,1:nFields)=(/1,4,2,3/)
00089          fdi(1,1:nFields)=(/1,3,4,2/)
00090          wdi(2,1:nFields)=(/1,4,3,2/)
00091          fdi(2,1:nFields)=(/1,4,3,2/)
00092       ELSEIF (nDim == 3) THEN
00093          nSteps=6
00094          ALLOCATE(dim_orders(nDim,nSteps))
00095          dim_orders=reshape((/1,2,3, 1,3,2, 2,1,3, 2,3,1, 3,1,2, 3,2,1/),(/nDim,nSteps/))
00096 
00097          nFields=5
00098          ALLOCATE(wdi(nDim,nFields+NrTracerVars), fdi(nDim,nFields+NrTracerVars))
00099          wdi(1,1:nFields)=(/1,5,2,3,4/)
00100          fdi(1,1:nFields)=(/1,3,4,5,2/)
00101          wdi(2,1:nFields)=(/1,5,3,4,2/)
00102          fdi(2,1:nFields)=(/1,5,3,4,2/)
00103          wdi(3,1:nFields)=(/1,5,4,2,3/)
00104          fdi(3,1:nFields)=(/1,4,5,3,2/)
00105       END IF
00106       IF (NrTracerVars > 0) THEN
00107          DO d=1,nDim
00108             wdi(d,nFields+1:nFields+NrTracerVars)=(/(nFields+i, i=1, NrTracerVars)/)
00109             fdi(d,nFields+1:nFields+NrTracerVars)=(/(nFields+i, i=1, NrTracerVars)/)
00110          END DO
00111          nFields=nFields+NrTracerVars
00112       END IF
00113 
00114       write(*,*) 'nFields=', nFields-NrTracerVars
00115       write(*,*) 'NrTracerVars=', NrTracerVars
00116       write(*,*) 'wdi=', wdi(1,:)
00117       write(*,*) 'fdi=', fdi(1,:)
00118       write(*,*) 'mhigh_1=', m_high+1
00119       ! stop
00120    END SUBROUTINE MUSCLInit
00121 
00122 SUBROUTINE MUSCLAdvance(Info, dom_range, lComplete)
00123   TYPE (InfoDef) :: Info
00124   LOGICAL :: lComplete
00125   INTEGER :: level
00126   INTEGER :: mx  !number of cells in data to update
00127   REAL(KIND=qPREC) :: dt, hdt !time and half time to advance for
00128   INTEGER, DIMENSION(3,2) :: dom_range, ghost_range
00129   INTEGER, DIMENSION(3,2) :: row, grow
00130   INTEGER, DIMENSION(3,2) :: ms 
00131   INTEGER :: d, i, j, k, m
00132   INTEGER ::  dir
00133   INTEGER , DIMENSION(2) :: pdirs
00134   REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: data !array that goes from 2-hyperbolic_mbc to mx+hyperbolic_mbc that needs to be updated from 1:mx
00135   REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: fluxes !array that has fluxes from 1:mx+1  
00136   REAL(KIND=qPREC) :: dtdx
00137 
00138   ghost_range=dom_range
00139 !  write(*,*) 'dom_range=', dom_range
00140   ghost_range(1:nDim,:)=ghost_range(1:nDim,:)+spread(hyperbolic_mbc*(/-1,1/),1,nDim)
00141 !  write(*,*) spread(hyperbolic_mbc*(/-1,1/),nDim,1)
00142 !  write(*,*) 'ghost range=', ghost_range
00143   dt=levels(Info%level)%dt
00144   hdt=.5d0*dt
00145 
00146   dtdx=dt/levels(Info%level)%dx
00147 
00148   CALL BeforeStep(Info)
00149   CALL Src(Info, ghost_range, levels(Info%level)%tnow, hdt)
00150 
00151   DO d=1,nDim
00152 
00153      ! direction of update
00154      dir=dim_orders(d,modulo(levels(Info%level)%CurrentLevelStep-1, nSteps)+1)
00155 
00156 !     write(*,*) 'doing pass for dimension ',dir, levels(Info%level)%CurrentLevelStep
00157 
00158      ! normal directions
00159      pdirs=modulo(dir+(/0,1/),3)+1
00160 
00161      ! grow is range of Info%q to map onto data including ghost zones
00162      grow(dir,:)=ghost_range(dir,:)
00163 
00164      ! row is range of Info%q to update
00165      row(dir,:)=dom_range(dir,:)
00166 
00167      ! ms is bounds for returned fluxes needed to update Info%q
00168      ms(dir,:)=dom_range(dir,:)+(/0,1/)
00169 
00170      ! mx is number of cells to update
00171      mx=row(dir,2)-row(dir,1)+1
00172 
00173      ALLOCATE(data(1-hyperbolic_mbc:mx+hyperbolic_mbc, nFields), fluxes(1:mx+1, nFields))
00174 
00175      DO j=ghost_range(pdirs(1),1), ghost_range(pdirs(1),2)
00176         !update bounds in transverse directions to update
00177         grow(pdirs(1),:)=j
00178         row(pdirs(1),:)=j
00179         mS(pdirs(1),:)=j
00180 
00181 !        write(*,*) 'j=', j
00182         DO k=ghost_range(pdirs(2),1), ghost_range(pdirs(2),2)
00183            !update bounds in transverse directions to update
00184            grow(pdirs(2),:)=k
00185            row(pdirs(2),:)=k
00186            mS(pdirs(2),:)=k
00187 
00188            !remap bounds of ghost zones onto data array and map fields in q to correct order (rho, E, v_perp, v_parallel... using wdi(:,:)
00189            data(1-hyperbolic_mbc:mx+hyperbolic_mbc,1:nFields)= &
00190                 reshape(Info%q(grow(1,1):grow(1,2), grow(2,1):grow(2,2), grow(3,1):grow(3,2), wDi(dir,:)),(/mx+2*hyperbolic_mbc,nFields/))
00191 
00192 
00193            !calculate fluxes using data
00194            CALL calc_fluxes(data, mx, levels(Info%level)%dt, levels(Info%level)%dx, fluxes, maxspeed(Info%level))
00195 
00196 !           FORALL(i=lbound(fluxes,1):ubound(fluxes, 1))
00197               !remap fluxes back to order of fields in q
00198               fluxes(:,:)=fluxes(:,fdi(dir,:))*dtdx
00199 !           end forall
00200 
00201            !pass fluxes to storefixupfluxes in case they need to be synchronized or restricted...
00202            CALL storefixupfluxes(Info,mS,dir,reshape(fluxes,(/ms(:,2)-ms(:,1)+1,nFields/)), (/(m,m=1,nFields)/))           
00203 
00204            !Update q with those fluxes
00205            Info%q(row(1,1):row(1,2),row(2,1):row(2,2),row(3,1):row(3,2),1:nFields)= &
00206                 Info%q(row(1,1):row(1,2),row(2,1):row(2,2),row(3,1):row(3,2),1:nFields)+ &
00207                 reshape((fluxes(1:mx,1:nFields)-fluxes(2:mx+1,1:nFields)), (/row(:,2)-row(:,1)+1,nFields/))
00208         END DO
00209      END DO
00210 
00211   !   DO i=1, Info%mX(dir)
00212   !      write(*,'(I,10E20.5)') i, Info%q(i,1,1,:)
00213   !   END DO
00214      !Now we can shrink the ghost_range in the direction we just updated
00215      ghost_range(dir,:)=ghost_range(dir,:)+hyperbolic_mbc*(/1,-1/)
00216 
00217      !free up work arrays for this direction
00218      deallocate(data, fluxes)
00219   END DO
00220 
00221   CALL Src(Info, dom_range, levels(Info%level)%tnow+hdt, hdt)
00222   
00223   CALL AfterStep(Info)
00224 
00225 END SUBROUTINE MUSCLAdvance
00226 
00227 subroutine calc_fluxes(data, mx, dt, dx, fluxes, maxspeed)
00228   
00229   REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: data, fluxes
00230   REAL(KIND=qPREC), POINTER, DIMENSION(:,:) :: w, dw
00231   REAL(KIND=qPREC), POINTER, DIMENSION(:) :: wl, wr, fl, fr, df, qr, ql
00232   INTEGER :: mx  !number of cells in data to update
00233   INTEGER :: i, mbc, j
00234   REAL(KIND=qPREC) :: dt, dx, dtdx, hdt  !time and half time to advance for
00235   REAL(KIND=qPREC), DIMENSION(:), ALLOCATABLE :: f
00236   REAL(KIND=qPREC) :: maxspeed, dright, dleft, dcenter
00237   !order of fields in data is rho, E, px, [py, [pz]]
00238   mbc=hyperbolic_mbc
00239   ALLOCATE(w(1-mbc:mx+mbc, nFields))
00240   ALLOCATE(dw(1-mbc+1:mx+mbc-1, nFields))
00241   !ALLOCATE(wl(1-mbc:mx+mbc))
00242   !ALLOCATE(wr(1-mbc+1:mx+mbc+1))
00243   ALLOCATE(wl(nFields))
00244   ALLOCATE(wr(nFields))
00245   ALLOCATE(f(nFields))
00246   ALLOCATE(fr(nFields))
00247   ALLOCATE(fl(nFields))
00248   ALLOCATE(df(nFields))
00249   ALLOCATE(qr(nFields))
00250   ALLOCATE(ql(nFields))
00251 
00252   dtdx=dt/dx
00253 
00254   DO i=1-mbc,mx+mbc
00255      CALL cons_to_prim_MUSCL(data(i,:), w(i,:))
00256   END DO
00257 
00258   DO i=1-mbc+1,mx+mbc-1
00259      IF (lUseLimiter) THEN
00260         DO j=1,nFields
00261            dright=w(i+1,j)-w(i,j)
00262            dleft=w(i,j)-w(i-1,j)
00263            dcenter=.5d0*(dright+dleft)
00264            IF (SIGN(1d0, dleft) == SIGN(1d0, dright)) THEN
00265               dw(i,j)=sign(min(2d0*abs(dleft), 2d0*abs(dright), abs(dcenter)), dcenter)
00266            ELSE
00267               dw(i,j)=0d0           
00268            END IF
00269         END DO
00270      ELSE
00271         dw(i,:)=.5d0*(w(i+1,:)-w(i-1,:))
00272      END IF
00273 !     dw(i,:)=0d0
00274 !     write(*,'(I,10E20.5)') i, dw(i,:)
00275   END DO
00276 
00277   DO i=1,mx+1
00278      wl = w(i-1,:)+.5d0*dw(i-1,:)
00279      wr = w(i,:)-.5d0*dw(i,:)
00280      maxspeed=max(maxspeed, calc_flux(wl, wr, f))
00281 !     write(*,'(A,I)') 'i = ', i
00282 !     write(*,'(A,10E20.5)') 'w(i-1) = ',w(i-1,:) 
00283 !     write(*,'(A,10E20.5)') 'dw(i-1) = ',dw(i-1,:) 
00284 !     write(*,'(A,10E20.5)') 'wl = ', wl
00285 !     write(*,'(A,10E20.5)') 'w(i) = ',w(i,:) 
00286 !     write(*,'(A,10E20.5)') 'dw(i) = ',dw(i,:) 
00287 !     write(*,'(A,10E20.5)') 'wr = ', wr
00288 !     write(*,'(A,10E20.5)') 'f = ', f
00289      fluxes(i,:)=f
00290   END DO
00291 
00292   DEALLOCATE(w, dw, wl, wr, f, fr, fl, df, qr, ql)
00293 end subroutine calc_fluxes 
00294 
00295 
00296 subroutine cons_to_prim_MUSCL(q,w)
00297   REAL(KIND=qPREC),DIMENSION(:) :: q,w
00298   w = q
00299   w(2)=gamma1*(q(2)-half*sum(q(3:m_high+1)**2)/q(1))
00300   w(3:m_high+1)=q(3:m_high+1)/q(1)
00301   !w(nDim+2+1:nFields)=q(nDim+2+1:nFields)
00302    
00303 END subroutine cons_to_prim_MUSCL
00304 
00305 subroutine prim_to_cons_MUSCL(w,q)
00306   REAL(KIND=qPREC),DIMENSION(:) :: w,q
00307   q = w
00308   q(2)=gamma7*w(2)+half*sum(w(3:m_high+1)**2)*w(1)
00309   q(3:m_high+1)=w(3:m_high+1)*w(1)
00310   !w(nDim+2+1:nFields)=q(nDim+2+1:nFields)
00311    
00312 END subroutine prim_to_cons_MUSCL
00313 
00314 END MODULE MUSCLScheme
 All Classes Files Functions Variables