Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:41:03 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_OSC_MUL_Y(iy,hh,mask,ohat,scal)
                0004 C     |================================================================|
                0005 C     | OSC_MUL_Y: evaluate WENO oscillation weights in Y.             |
                0006 C     |================================================================|
                0007 
                0008           implicit none
                0009 
                0010 C     =============================================== global variables
                0011 #         include "SIZE.h"
                0012 
                0013           integer iy,hh
                0014           _RL mask(1-OLy:sNy+OLy)
                0015           _RL ohat(1:2,
                0016      &             1-OLy:sNy+OLy)
                0017           _RL scal(1:2)
                0018 
                0019           integer ii
                0020           _RL dels,dfs1,dfs2
                0021           _RL osum,zero,mval
                0022           _RL oval,omin,omax
                0023 
                0024 C     =============================== calc. WENO oscillation weighting
                0025 c         omin = +huge(+1. _d 0)
                0026 c         omax = -huge(+1. _d 0)
                0027           omin = +1. _d 99
                0028           omax = -1. _d 99
                0029 
                0030           zero = 1. _d -20
                0031           mval = 1. _d + 0
                0032 
                0033           do  ii = iy-hh, iy+hh
                0034 
                0035 C     =============================== calc. derivatives centred on II.
                0036               dels = (ii - iy) * 2. _d 0
                0037 
                0038               dfs1 = ohat(1,ii)
                0039               dfs2 = ohat(2,ii)
                0040 
                0041               dfs1 = dfs1 + dfs2 * dels
                0042 
                0043 C     =============================== oscl. = NORM(H^N * D^N/DY^N(F)).
                0044               oval = (2. _d 0 * dfs1)**2
                0045      &             + (4. _d 0 * dfs2)**2
                0046 
                0047               if (oval.lt.omin) omin = oval
                0048               if (oval.gt.omax) omax = oval
                0049 
                0050 C     =============================== any mask across oscil. stencil
                0051               mval = mval * mask(ii)
                0052 
                0053           end do
                0054 
                0055           if (mval .gt. 0. _d 0) then
                0056 
                0057 C     =============================== calc. WENO-style profile weights
                0058               scal(1) = 1. _d 5
                0059      &             / (omax + zero)**3
                0060               scal(2) = 1. _d 0
                0061      &             / (omin + zero)**3
                0062 
                0063               osum = scal(1) + scal(2)
                0064               scal(1) = scal(1) / osum
                0065               scal(2) = scal(2) / osum
                0066 
                0067           else
                0068 
                0069 C     =============================== default to MONO. profile weights
                0070               scal(1) = 0. _d 0
                0071               scal(2) = 1. _d 0
                0072 
                0073           end if
                0074 
                0075           return
                0076 
                0077 c     end subroutine GAD_OSC_MUL_Y
                0078       end