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