Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:41:04 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_PPM_HAT_Y(bi,bj,kk,ix,
                0004      I           method,mask,fbar,edge,
                0005      I           ohat,fhat,myThid)
                0006 C     |================================================================|
                0007 C     | PPM_HAT_Y: reconstruct grid-cell PPM polynomials.              |
                0008 C     |================================================================|
                0009 
                0010           implicit none
                0011 
                0012 C     =============================================== global variables
                0013 #         include "SIZE.h"
                0014 #         include "GRID.h"
                0015 #         include "GAD.h"
                0016 
                0017 C     ====================================================== arguments
                0018 C       bi,bj     :: tile indexing.
                0019 C       kk        :: r-indexing.
                0020 C       ix        :: x-indexing.
                0021 C       method    :: advection scheme.
                0022 C       mask      :: row of cell-wise masking values.
                0023 C       fbar      :: row of cell-wise values.
                0024 C       edge      :: row of edge-wise values.
                0025 C       ohat      :: row of oscl. coeff.
                0026 C                    - OHAT(1,:) = D^1F/DS^1.
                0027 C                    - OHAT(2,:) = D^2F/DS^2.
                0028 C       fhat      :: row of poly. coeff.
                0029 C                    - FHAT(:,I) = PPM coeff.
                0030 C       myThid    :: thread number.
                0031 C     ================================================================
                0032           integer bi,bj,kk,ix
                0033           integer method
                0034           _RL mask(1-OLy:sNy+OLy)
                0035           _RL fbar(1-OLy:sNy+OLy)
                0036           _RL edge(1-OLy:sNy+OLy)
                0037           _RL ohat(1:2,
                0038      &             1-OLy:sNy+OLy)
                0039           _RL fhat(1:3,
                0040      &             1-OLy:sNy+OLy)
                0041           integer myThid
                0042 
                0043 C     ====================================================== variables
                0044 C       ii,iy     :: local y-indexing.
                0045 C       ff00      :: centre-biased cell mean.
                0046 C       ffll,ffrr :: left-, and right-biased cell values.
                0047 C       fell,ferr :: left-, and right-biased edge values.
                0048 C       dfds      :: linear slope estimates.
                0049 C       uhat      :: "NULL" limited poly. coeff.
                0050 C       lhat      :: "MONO" limited poly. coeff.
                0051 C       scal      :: "WENO" weights.
                0052 C       fmag,fdel :: local perturbation indicators.
                0053 C     ================================================================
                0054           integer ii,iy
                0055           _RL ff00,ffll,ffrr
                0056           _RL fell,ferr
                0057           _RL dfds(-1:+1)
                0058           _RL uhat(+1:+5)
                0059           _RL lhat(+1:+5)
                0060           _RL scal(+1:+2)
                0061           _RL fmag,fdel
                0062           integer  mono
                0063 
                0064 C     ==================== reconstruct coeff. for grid-cell poynomials
                0065           do  iy = 1-OLy+2, sNy+OLy-2
                0066 
                0067 C     =============================== assemble cell mean + edge values
                0068               ff00 = fbar(iy+0)
                0069               ffll = ff00
                0070      &             + mask(iy-1)*(fbar(iy-1)-ff00)
                0071               ffrr = ff00
                0072      &             + mask(iy+1)*(fbar(iy+1)-ff00)
                0073 
                0074               fell = edge(iy-0)
                0075               ferr = edge(iy+1)
                0076 
                0077 c             select case(method)
                0078 c             case(ENUM_PPM_NULL_LIMIT)
                0079               if ( method.eq.ENUM_PPM_NULL_LIMIT ) then
                0080 C     =============================== "NULL" limited grid-cell profile
                0081               CALL GAD_PPM_FUN_NULL ( ff00,fell,ferr,
                0082      &             lhat,mono)
                0083 
                0084 c             case(ENUM_PPM_MONO_LIMIT)
                0085               elseif ( method.eq.ENUM_PPM_MONO_LIMIT ) then
                0086 C     =============================== "MONO" limited grid-cell profile
                0087               CALL GAD_PLM_FUN_U(ffll,ff00,ffrr,dfds)
                0088 
                0089               CALL GAD_PPM_FUN_MONO ( ff00,ffll,ffrr,
                0090      &             fell,ferr,dfds,lhat,mono)
                0091 
                0092 c             case(ENUM_PPM_WENO_LIMIT)
                0093               elseif ( method.eq.ENUM_PPM_WENO_LIMIT ) then
                0094 C     =============================== "WENO" limited grid-cell profile
                0095               CALL GAD_PLM_FUN_U(ffll,ff00,ffrr,dfds)
                0096 
                0097               CALL GAD_PPM_FUN_NULL ( ff00,fell,ferr,
                0098      &             uhat,mono)
                0099 
                0100               CALL GAD_PPM_FUN_MONO ( ff00,ffll,ffrr,
                0101      &             fell,ferr,dfds,lhat,mono)
                0102 
                0103               if ( mono .gt. 0) then
                0104 
                0105 C     =============================== only apply WENO if it is worth it
                0106               fdel = abs(ffrr-ff00)+abs(ff00-ffll)
                0107               fmag = abs(ffll)+abs(ff00)+abs(ffrr)
                0108 
                0109               if (fdel .gt. 1. _d -6 * fmag) then
                0110 
                0111 C     =============================== calc. WENO oscillation weighting
                0112                   CALL GAD_OSC_MUL_Y(iy,+2,mask,
                0113      &                 ohat,scal)
                0114 
                0115                   do  ii = +1, +3
                0116 C     =============================== blend limited/un-limited profile
                0117                   lhat(ii) = scal(1) * uhat(ii)
                0118      &                     + scal(2) * lhat(ii)
                0119                   end do
                0120 
                0121               end if
                0122 
                0123               end if
                0124 
                0125 c             end select
                0126               endif
                0127 
                0128               do  ii = +1, +3
                0129 C     =============================== copy polynomial onto output data
                0130                   fhat(ii,iy) = lhat(ii)
                0131               end do
                0132 
                0133           end do
                0134 
                0135           return
                0136 
                0137 c     end subroutine GAD_PPM_HAT_Y
                0138       end