Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:41:05 UTC

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
8e4c181d69 Jean*0001 #     include "GAD_OPTIONS.h"
                0002 
                0003       SUBROUTINE GAD_PQM_ADV_Y(meth,bi,bj,kk,
                0004      I           calc_CFL,delT,vvel,vfac,fbar,
                0005      O           flux,myThid )
                0006 C     |================================================================|
                0007 C     | PQM_ADV_Y: evaluate grid-cell advective flux in Y.             |
                0008 C     | Lagrangian-type Piecewise Quartic Method (PQM).                |
                0009 C     |================================================================|
                0010 
                0011           implicit none
                0012 
                0013 C     =============================================== global variables
                0014 #         include "SIZE.h"
                0015 #         include "GRID.h"
                0016 #         include "GAD.h"
                0017 
                0018 C     ================================================================
                0019 C       meth     :: advection method.
                0020 C       bi,bj    :: tile indexing.
                0021 C       kk       :: r-index.
                0022 C       calc_CFL :: TRUE to calc. CFL from vel.
                0023 C       delT     :: time-step.
                0024 C       vvel     :: vel.-comp in y-direction.
                0025 C       vfac     :: vel.-flux in y-direction.
                0026 C       fbar     :: grid-cell values.
                0027 C       flux     :: adv.-flux in y-direction.
                0028 C       myThid   :: thread number.
                0029 C     ================================================================
                0030           integer meth
                0031           integer bi,bj,kk
                0032           logical calc_CFL
                0033           _RL delT
                0034           _RL vvel(1-OLx:sNx+OLx,
                0035      &             1-OLy:sNy+OLy)
                0036           _RL vfac(1-OLx:sNx+OLx,
                0037      &             1-OLy:sNy+OLy)
                0038           _RL fbar(1-OLx:sNx+OLx,
                0039      &             1-OLy:sNy+OLy)
                0040           _RL flux(1-OLx:sNx+OLx,
                0041      &             1-OLy:sNy+OLy)
                0042           integer myThid
                0043 
                0044 C     ================================================================
                0045 C       ix,iy,ir :: grid indexing.
                0046 C       floc     :: row of grid-cell values.
                0047 C       mloc     :: row of grid-cell mask values.
                0048 C       fhat     :: row of poly. coeff.
                0049 C                    - FHAT(:,I) = PQM coeff.
                0050 C       edge     :: row of edge-wise values/slopes.
                0051 C                    - EDGE(1,:) = VALUE.
                0052 C                    - EDGE(2,:) = DF/DY.
                0053 C       ohat     :: row of oscl. coeff.
                0054 C                    - OHAT(1,:) = D^1F/DS^1.
                0055 C                    - OHAT(2,:) = D^2F/DS^2.
                0056 C     ================================================================
                0057           integer ix,iy
                0058           _RL mloc(1-OLy:sNy+OLy)
                0059           _RL floc(1-OLy:sNy+OLy)
                0060           _RL fhat(1:5,
                0061      &             1-OLy:sNy+OLy)
                0062           _RL edge(1:2,
                0063      &             1-OLy:sNy+OLy)
                0064           _RL ohat(1:2,
                0065      &             1-OLy:sNy+OLy)
                0066           _RL vsum
                0067 
                0068           do ix = 1-OLx+0, sNx+OLx-0
                0069 C     ==================== zero stencil "ghost" cells along boundaries
                0070               flux(ix, +1-OLy+0) = 0. _d 0
                0071               flux(ix, +1-OLy+1) = 0. _d 0
                0072               flux(ix, +1-OLy+2) = 0. _d 0
                0073               flux(ix, +1-OLy+3) = 0. _d 0
                0074               flux(ix,sNy+OLy-0) = 0. _d 0
                0075               flux(ix,sNy+OLy-1) = 0. _d 0
                0076               flux(ix,sNy+OLy-2) = 0. _d 0
                0077           end do
                0078 
                0079 C     ================================================================
                0080 C       (1): copy a single row of data onto contiguous storage, treat
                0081 C            as a set of one-dimensional problems.
                0082 C       (2): calc. "oscillation-indicators" for each grid-cell if ad-
                0083 C            vection scheme is WENO-class.
                0084 C       (3): calc. edge-centred values/slopes by high-order interpol-
                0085 C            ation.
                0086 C       (4): calc. cell-centred polynomial profiles with appropriate
                0087 C            slope-limiting.
                0088 C       (5): calc. fluxes using a local, semi-lagrangian integration.
                0089 C     ================================================================
                0090 
                0091           do ix = 1-OLx+0, sNx+OLx-0
                0092 
                0093           vsum = 0.0 _d 0
                0094           do iy = 1-OLy+0, sNy+OLy-0
                0095 C     ================================== quick break on zero transport
                0096               vsum = vsum
                0097      &             + abs(vfac(ix,iy))
                0098 C     ================================== make local unit-stride copies
                0099               floc(iy) = fbar (ix,iy)
                0100               mloc(iy) =
                0101      &          maskC(ix,iy,kk,bi,bj)
                0102           end do
                0103 
                0104           if (vsum .gt. 0. _d 0) then
                0105 
                0106 C     ==================== reconstruct derivatives for WENO indicators
                0107           if (meth.eq.ENUM_PQM_WENO_LIMIT) then
                0108           CALL GAD_OSC_HAT_Y(bi,bj,kk,ix,
                0109      &                   mloc,floc,
                0110      &                   ohat,myThid)
                0111           end if
                0112 
                0113 C     ==================== reconstruct 5th--order accurate edge values
                0114           CALL GAD_PQM_P5E_Y(bi,bj,kk,ix,
                0115      &                   mloc,floc,
                0116      &                   edge,myThid)
                0117 
                0118 C     ==================== reconstruct coeff. for grid-cell poynomials
                0119           CALL GAD_PQM_HAT_Y(bi,bj,kk,ix,
                0120      &                   meth,
                0121      &                   mloc,floc,
                0122      &                   edge,ohat,
                0123      &                   fhat,myThid)
                0124 
                0125 C     ==================== evaluate integral fluxes on grid-cell edges
                0126           CALL GAD_PQM_FLX_Y(bi,bj,kk,ix,
                0127      &                   calc_CFL,
                0128      &                   delT,vvel,
                0129      &                   vfac,fhat,
                0130      &                   flux,myThid)
                0131 
                0132           else
                0133 
                0134           do iy = 1-OLy+4, sNy+OLy-3
                0135 C     ================================== "null" flux on zero transport
                0136               flux(ix,iy) = 0.0 _d 0
                0137           end do
                0138 
                0139           end if
                0140 
                0141           end do
                0142 
                0143           return
                0144 
                0145 c     end subroutine GAD_PQM_ADV_Y
                0146       end