• VTF
  • FSI
  • AMROC
  • SFC
  • Motion
  • STLIB
  • Main Page
  • src/generic/mesh.f90

    
    MODULE mesh
      SAVE
    
      ! ----  the declarationf for 
      ! ----  the upper/right (r) and lower/left (l) boundaries
      ! ----  of our (up to) 3D patch including the ghost-cells
    
      ! ----  The grid spacing is (dx,dy,dz)
      ! ----  The cells are centered at (xc,yc,zc)
      
      
      DOUBLE PRECISION:: xl, xr
      DOUBLE PRECISION:: yl, yr
      DOUBLE PRECISION:: zl, zr
      DOUBLE PRECISION:: dx, dy, dz
    
      ! ----  the declarationf for 
      ! ----  the upper/right (r) and lower/left (l) boundaries
      ! ----  of the entire (up to) 3D domain excluding the ghost-cells
    
      DOUBLE PRECISION:: xlg, xrg
      DOUBLE PRECISION:: ylg, yrg
      DOUBLE PRECISION:: zlg, zrg
    
    END MODULE mesh
    
    
    MODULE array_bounds
      SAVE
    
      ! ----  The array indicies of the upper(hi) and lower (lo)
      ! ----  bounds of our domain, including ghost-cells
    
    
      ! ----  The array address of the physical domain 
      ! ----  is (1:nx) x (1:ny) x (1:nz)
    
    
      ! ----  Most of these are set in Generic_SetArrayBounds.f90
    
      INTEGER:: nx,ny,nz
      INTEGER:: mx,my,mz
    
      ! ----  Bounds including ghost-cells for arrays
    
      INTEGER:: ixlo
      INTEGER:: ixhi
    
      INTEGER:: iylo
      INTEGER:: iyhi
     
      INTEGER:: izlo
      INTEGER:: izhi
    
      ! ----  Generic array bounds used when doing x,y or z sweeps.
      INTEGER:: inlo
      INTEGER:: inhi
    
      ! ----  Generic array bounds used when doing x,y or sweeps.
      ! ----  across derivatives
      INTEGER:: idnlo
      INTEGER:: idnhi
    
    END MODULE array_bounds
    
    
    
    MODULE method_parms
      SAVE
    
      ! ----  The parameters that relate
      ! ----  the number of ghost cells to the size of a 
      ! ----  given ENO stencil.
       
      ! ----  The parameter nvar is the number of 
      ! ----  elements in the state vector - 
      ! ----  This number should be >= 6 ( for the 3D problem)
      ! ----  5 for the phsycial dependant variables, 1 for a
      ! ----  passive scalar used to determine gamma
      ! ----  more for other passive scalars
    
      INTEGER :: dim
      INTEGER :: ncomps
      INTEGER :: nvars
      INTEGER :: nscal
      INTEGER :: enoOrder 
      INTEGER :: nghost 
      
      ! ----  Flag set for using standard scheme 
      ! ----  or optimized ( set to 1 for optimized)
      INTEGER :: order
      INTEGER :: optimized
      INTEGER :: use_carbfix
      INTEGER :: use_dcflag
      INTEGER :: useExOutput 
    
      ! ----  Determins the width of the stencil to use.
      ! ----  enoOder = (stencil - 1) /2 is the order of 
      ! ----  accuracy of a canidate stencil inside WENO
      INTEGER :: stencil
    
      ! ----  upper and lower bounds for center difference interp
      INTEGER :: upb,lob
      
    
      ! ----  Selects the correct different solver
      ! ----  method = 0 is WENO
      ! ----  method = 1 is WENO-TCD
      ! ----  method = 2 is 'New Hybrid'
      INTEGER :: method = 0
    
      ! ---- if not using Time Refinement set to 1.
      INTEGER :: noTimeRefine = 0
      
      ! ----  Select the resolved scale viscous terms
      ! ----  useViscous = 0 no viscous terms
      ! ----  useViscous = 1 use viscous terms
      ! ----       ?      this requires more ghost cells
      INTEGER :: useViscous = 0
    
    
      ! ----  Determines if LES is being used in 3D
      ! ----  useLES = 0 no LES model
      ! ----  useLES = 1 use the LES model 
      ! ----             this requires more ghost cells
      INTEGER :: useLES = 0
      
      ! ---- Determine if the source term is to be used
      ! ---- useSource = 0 no source term
      ! ---- useSource = 1 source term - defined as weno_source 
      !                    in local problem directory
      INTEGER :: useSource = 0
      
      ! ---- Are periodic boundary conditions used 
      ! ---- in the x,y,z direction
      INTEGER :: xper = 0
      INTEGER :: yper = 0
      INTEGER :: zper = 0
    
      ! these flags control the characteristic boundary condition type
      INTEGER :: cbc_direction(6)
      DATA cbc_direction /0, 0, 0, 0, 0, 0/
      INTEGER :: cbc_ixlow, cbc_ixup
      INTEGER :: cbc_iylow, cbc_iyup
      INTEGER :: cbc_izlow, cbc_izup
    
      ! these flags control the boundary stencil
      INTEGER :: bc_ixlow, bc_ixup
      INTEGER :: bc_iylow, bc_iyup
      INTEGER :: bc_izlow, bc_izup
    
      ! filtering variables
      DOUBLE PRECISION :: alpha_eta1, alpha_eta2
    
    END MODULE method_parms
    
    subroutine cles_getiparam(name, value, ierr)
      
      use method_parms
    
      implicit none
      
      include 'cles.i'
    
      character(LEN=*) :: name
      integer :: value, ierr
      
      character(LEN=12) :: cname
    
      cname = name
      ierr = CLES_RETURN_OK
    
      select case (cname)
      case ('nvars')
         value = nvars
      case ('ncomps')
         value = ncomps
      case ('dim')
         value = dim
      case ('nscal')
         value = nscal
      case ('stencil')
         value = stencil
      case ('method')
         value = method
      case ('notimerefine')
         value = noTimeRefine
      case ('order')
         value = order
      case ('nghost')
         value = nghost
      case ('optimized')
         value = optimized
      case ('use_carbfix')
         value = use_carbfix
      case ('enoorder')
         value = enoOrder
      case ('useles')
         value = useLES
      case ('useviscous')
         value = useViscous
      case ('usesource')
         value = useSource
      case default
         ierr = CLES_RETURN_ERROR
      end select
    
    end subroutine cles_getiparam
    
    MODULE Interp_coeffs
      
      USE method_parms
      
      SAVE
      
      ! ----  The declarations for the coefficients used 
      ! ----  in calculating flux interpolants by WENO
      
      ! ---- make these too big, but static
      DOUBLE PRECISION :: aw(0:10,0:10-1)
      DOUBLE PRECISION :: cw(0:10)
      DOUBLE PRECISION :: dw(0:10,10-1,0:10-1)
      
      DOUBLE PRECISION :: ac, cc, dc, b2, b3
    
      ! ----  used in the Center Difference (Finite difference)
      INTEGER :: tcd_bndry_width = 2
      INTEGER :: tcd_bndry_points = 4
      DOUBLE PRECISION :: tcd_center(3)
      DOUBLE PRECISION :: tcd_flux(3)
      ! dimensions of tcd_bndry are max number of points (4) x max number of stencils (2) 
      DOUBLE PRECISION :: tcd_bndry(4,2)
    
    END MODULE Interp_coeffs
    
    module cles_interfaces
    
      interface 
         subroutine cles_xLocation(i, x)
           integer i
           double precision x
         end subroutine cles_xLocation
      end interface
    
      interface 
         subroutine cles_yLocation(i, x)
           integer i
           double precision x
         end subroutine cles_yLocation
      end interface
    
      interface 
         subroutine cles_zLocation(i, x)
           integer i
           double precision x
         end subroutine cles_zLocation
      end interface
    
      interface
         subroutine cles_eqstate(q,ncomps,qt,nvars,n,useLES)
           integer  n, ncomps, nvars, useLES
           DOUBLE PRECISION  q(ncomps,n), qt(nvars,n)
         end subroutine cles_eqstate
      end interface
    
      interface 
         subroutine cles_roe(ux0, ux1, ncomps, vx0, vx1, nvars, &
              gamma, Rr, mu, nscal)
           integer ncomps, nvars, nscal
           double precision ux0(ncomps), ux1(ncomps)
           double precision vx0(nvars), vx1(nvars)
           double precision gamma, Rr, mu(*)
         end subroutine cles_roe
      end interface
    
      interface 
         subroutine SetUpLES(version_)
           INTEGER, INTENT(IN) :: version_
         end subroutine SetUpLES
      end interface
    
      interface
         subroutine InitializeLES(ux,vx)
           use mesh
           use array_bounds
           use method_parms
    
           DOUBLE PRECISION, INTENT(INOUT) :: ux(ncomps,ixlo:ixhi,iylo:iyhi,izlo:izhi)
           DOUBLE PRECISION, INTENT(IN) :: vx(nvars,ixlo:ixhi,iylo:iyhi,izlo:izhi)
         end subroutine InitializeLES
      end interface
    
      interface 
         subroutine AddSgsFlux(fxi,ux,vx,direction)
           use mesh
           use array_bounds
           use method_parms
    
           DOUBLE PRECISION, INTENT(IN) :: ux(ncomps,ixlo:ixhi,iylo:iyhi,izlo:izhi)
           DOUBLE PRECISION, INTENT(IN) :: vx(nvars,ixlo:ixhi,iylo:iyhi,izlo:izhi)
           DOUBLE PRECISION, INTENT(INOUT)::fxi(ncomps,ixlo:ixhi,iylo:iyhi,izlo:izhi,3)
           INTEGER, INTENT(IN) :: direction
         end subroutine AddSgsFlux
      end interface
    
      interface 
         integer function cleslog_active(mod, loc)
           integer :: mod
           integer :: loc
         end function cleslog_active
      end interface
    
      interface 
         subroutine cleslog_log(slog)
           character(LEN=*) :: slog
         end subroutine cleslog_log
      end interface
    
      interface 
         subroutine cleslog_log_flush()
         end subroutine cleslog_log_flush
      end interface
    
    end module cles_interfaces
    

<