|
||||
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 UTC8e4c181d69 Jean*0001 # include "GAD_OPTIONS.h" 0002 0003 SUBROUTINE GAD_PPM_ADV_R(meth,bi,bj, 0004 I delT,velR,facR,fbar, 0005 O flux,myThid ) 0006 C |================================================================| 0007 C | PPM_ADV_R: evaluate grid-cell advective flux in R. | 0008 C | Lagrangian-type Piecewise Parabolic Method (PPM). | 0009 C |================================================================| 0010 0011 implicit none 0012 0013 C =============================================== global variables 0014 # include "SIZE.h" 0015 # include "GRID.h" 0016 # include "GAD.h" 0017 0018 C ================================================================ 0019 C meth :: advection method. 0020 C bi,bj :: tile indexing 0021 C delT :: level-wise time-steps. 0022 C velR :: vel. field in r-direction. 0023 C facR :: grid-areas in r-direction. 0024 C fbar :: grid-cell values. 0025 C flux :: trac.-flux in r-direction. 0026 C myThid :: thread number. 0027 C ================================================================ 0028 integer meth 0029 integer bi,bj 0030 _RL delT(1:Nr) 0031 _RL velR(1-OLx:sNx+OLx, 0032 & 1-OLy:sNy+OLy, 0033 & 1:Nr) 0034 _RL facR(1-OLx:sNx+OLx, 0035 & 1-OLy:sNy+OLy, 0036 & 1:Nr) 0037 _RL fbar(1-OLx:sNx+OLx, 0038 & 1-OLy:sNy+OLy, 0039 & 1:Nr) 0040 _RL flux(1-OLx:sNx+OLx, 0041 & 1-OLy:sNy+OLy, 0042 & 1:Nr) 0043 integer myThid 0044 0045 C ================================================================ 0046 C ix,iy,ir :: grid indexing. 0047 C floc :: col. of grid-cell values. 0048 C mloc :: col. of grid-cell mask values. 0049 C fhat :: col. of poly. coeff. 0050 C - FHAT(:,I) = PQM coeff. 0051 C edge :: col. of edge-wise values/slopes. 0052 C - EDGE(1,:) = VALUE. 0053 C - EDGE(2,:) = DF/DR. 0054 C ohat :: col. of oscl. coeff. 0055 C - OHAT(1,:) = D^1F/DS^1. 0056 C - OHAT(2,:) = D^2F/DS^2. 0057 C ================================================================ 0058 integer ix,iy,ir 0059 _RL floc( 1-3:Nr+3) 0060 _RL mloc( 1-3:Nr+3) 0061 _RL fhat(1:3,1-0:Nr+0) 0062 _RL edge( 1-0:Nr+1) 0063 _RL ohat(1:2,1-3:Nr+3) 0064 _RL vsum 0065 0066 C ======================================= mask boundary conditions 0067 mloc( -2) = 0.0 _d 0 0068 mloc( -1) = 0.0 _d 0 0069 mloc( +0) = 0.0 _d 0 0070 mloc(Nr+1) = 0.0 _d 0 0071 mloc(Nr+2) = 0.0 _d 0 0072 mloc(Nr+3) = 0.0 _d 0 0073 0074 C ================================================================ 0075 C (1): copy a single row of data onto contiguous storage, treat 0076 C as a set of one-dimensional problems. 0077 C (2): calc. "oscillation-indicators" for each grid-cell if ad- 0078 C vection scheme is WENO-class. 0079 C (3): calc. edge-centred values/slopes by high-order interpol- 0080 C ation. 0081 C (4): calc. cell-centred polynomial profiles with appropriate 0082 C slope-limiting. 0083 C (5): calc. fluxes using a local, semi-lagrangian integration. 0084 C ================================================================ 0085 0086 do iy = 1-OLy+0, sNy+OLy-0 0087 do ix = 1-OLx+0, sNx+OLx-0 0088 C ======================================= no flux through surf. bc 0089 flux(ix,iy,+1) = 0.0 _d +0 0090 end do 0091 end do 0092 0093 C ==================== calculate transport for interior grid-cells 0094 do iy = 1-OLy+0, sNy+OLy-0 0095 do ix = 1-OLx+0, sNx+OLx-0 0096 0097 vsum = 0.0 _d 0 0098 do ir = 2, Nr 0099 C ================================== quick break on zero transport 0100 vsum = vsum 0101 & + abs(velR(ix,iy,ir)) 0102 end do 0103 0104 if (vsum .gt. 0. _d 0) then 0105 0106 do ir = 1, Nr 0107 C ================================== make local unit-stride copies 0108 floc(ir) = fbar(ix,iy,ir) 0109 mloc(ir) = 0110 & maskC(ix,iy,ir,bi,bj) 0111 end do 0112 0113 C ================================== make mask boundary conditions 0114 floc( -2) = floc(+1) 0115 floc( -1) = floc(+1) 0116 floc( +0) = floc(+1) 0117 floc(Nr+1) = floc(Nr) 0118 floc(Nr+2) = floc(Nr) 0119 floc(Nr+3) = floc(Nr) 0120 0121 C ==================== reconstruct derivatives for WENO indicators 0122 if (meth.eq.ENUM_PPM_WENO_LIMIT) then 0123 CALL GAD_OSC_HAT_R(bi,bj,ix,iy, 0124 & mloc,floc, 0125 & ohat,myThid) 0126 end if 0127 0128 C ==================== reconstruct 3rd--order accurate edge values 0129 CALL GAD_PPM_P3E_R(bi,bj,ix,iy, 0130 & mloc,floc, 0131 & edge,myThid) 0132 0133 C ==================== reconstruct coeff. for grid-cell poynomials 0134 CALL GAD_PPM_HAT_R(bi,bj,ix,iy, 0135 & meth, 0136 & mloc,floc, 0137 & edge,ohat, 0138 & fhat,myThid) 0139 0140 C ==================== evaluate integral fluxes on grid-cell edges 0141 CALL GAD_PPM_FLX_R(bi,bj,ix,iy, 0142 & delT,velR, 0143 & facR,fhat, 0144 & flux,myThid) 0145 0146 else 0147 0148 do ir = 2, Nr 0149 C ================================== "null" flux on zero transport 0150 flux(ix,iy,ir) = 0.0 _d +0 0151 end do 0152 0153 end if 0154 0155 end do 0156 end do 0157 0158 return 0159 0160 c end subroutine GAD_PPM_ADV_R 0161 end
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated from https://github.com/MITgcm/MITgcm by the 2.2.1-MITgcm-0.1 LXR engine. The LXR team |