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