Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:41:06 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_P5E_Y(bi,bj,kk,ix,
                0004      &           mask,fbar,edge,myThid)
                0005 C     |================================================================|
                0006 C     | PQM_P5E_Y: approximate edge values with degree-5 polynomials.  |
                0007 C     | Fixed grid-spacing variant in Y.                               |
                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           integer bi,bj,kk,ix
                0019           _RL mask(1-OLy:sNy+OLy)
                0020           _RL fbar(1-OLy:sNy+OLy)
                0021           _RL edge(1:2,
                0022      &             1-OLy:sNy+OLy)
                0023           integer myThid
                0024 
                0025 C     ====================================================== variables
                0026           integer iy
                0027           _RL mloc(-3:+2)
                0028           _RL floc(-3:+2)
                0029           _RL ftmp
                0030 
                0031 C     ==================== reconstruct 5th--order accurate edge values
                0032           do  iy = 1-OLy+3, sNy+OLy-2
                0033 
                0034 C     ================ mask local stencil: expand from centre outwards
                0035               mloc(-1) = mask(iy-1)
                0036               mloc(+0) = mask(iy+0)
                0037 
                0038               floc(-1) = fbar(iy+0)
                0039      &          + mloc(-1)*(fbar(iy-1)-fbar(iy+0))
                0040               floc(+0) = fbar(iy-1)
                0041      &          + mloc(+0)*(fbar(iy+0)-fbar(iy-1))
                0042 
                0043               mloc(-2) = mask(iy-2) * mloc(-1)
                0044               mloc(-3) = mask(iy-3) * mloc(-2)
                0045 
                0046               ftmp = 2. _d 0 * floc(-1) - floc(+0)
                0047               floc(-2) = ftmp
                0048      &          + mloc(-2)*(fbar(iy-2)-ftmp)
                0049               ftmp = 2. _d 0 * floc(-2) - floc(-1)
                0050               floc(-3) = ftmp
                0051      &          + mloc(-3)*(fbar(iy-3)-ftmp)
                0052 
                0053               mloc(+1) = mask(iy+1) * mloc(+0)
                0054               mloc(+2) = mask(iy+2) * mloc(+1)
                0055 
                0056               ftmp = 2. _d 0 * floc(+0) - floc(-1)
                0057               floc(+1) = ftmp
                0058      &          + mloc(+1)*(fbar(iy+1)-ftmp)
                0059               ftmp = 2. _d 0 * floc(+1) - floc(+0)
                0060               floc(+2) = ftmp
                0061      &          + mloc(+2)*(fbar(iy+2)-ftmp)
                0062 
                0063 C     ================ centred, 5th-order interpolation for edge value
                0064               edge(1,iy) =
                0065      &    +( 1. _d 0/60. _d 0)*(floc(-3)+floc(+2))
                0066      &    -( 8. _d 0/60. _d 0)*(floc(-2)+floc(+1))
                0067      &    +(37. _d 0/60. _d 0)*(floc(-1)+floc(+0))
                0068 
                0069               edge(2,iy) = (
                0070      &    -( 1. _d 0/90. _d 0)*(floc(-3)-floc(+2))
                0071      &    +( 5. _d 0/36. _d 0)*(floc(-2)-floc(+1))
                0072      &    -(49. _d 0/36. _d 0)*(floc(-1)-floc(+0))
                0073      &                ) * recip_dyC(ix,iy,bi,bj)
                0074 
                0075           end do
                0076 
                0077           return
                0078 
                0079 c     end subroutine GAD_PQM_P5E_Y
                0080       end