|
||||
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 UTC8e4c181d69 Jean*0001 # include "GAD_OPTIONS.h" 0002 0003 SUBROUTINE GAD_PQM_ADV_R(meth,bi,bj, 0004 I delT,wvel,wfac,fbar, 0005 O flux,myThid ) 0006 C |================================================================| 0007 C | PQM_ADV_R: evaluate grid-cell advective flux in R. | 0008 C | Lagrangian-type Piecewise Quartic Method (PQM). | 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 wvel :: vel.-comp in r-direction. 0023 C wfac :: vel.-flux in r-direction. 0024 C fbar :: grid-cell values. 0025 C flux :: adv.-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 wvel(1-OLx:sNx+OLx, 0032 & 1-OLy:sNy+OLy, 0033 & 1:Nr) 0034 _RL wfac(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:5,1-0:Nr+0) 0062 _RL edge(1:2,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 = +1, Nr 0099 C ================================== quick break on zero transport 0100 vsum = vsum 0101 & + abs(wfac(ix,iy,ir)) 0102 C ================================== make local unit-stride copies 0103 floc(ir) = fbar (ix,iy,ir) 0104 mloc(ir) = 0105 & maskC(ix,iy,ir,bi,bj) 0106 end do 0107 0108 if (vsum .gt. 0. _d 0) then 0109 0110 C ================================== make mask boundary conditions 0111 floc( -2) = floc(+1) 0112 floc( -1) = floc(+1) 0113 floc( +0) = floc(+1) 0114 floc(Nr+1) = floc(Nr) 0115 floc(Nr+2) = floc(Nr) 0116 floc(Nr+3) = floc(Nr) 0117 0118 C ==================== reconstruct derivatives for WENO indicators 0119 if (meth.eq.ENUM_PQM_WENO_LIMIT) then 0120 CALL GAD_OSC_HAT_R(bi,bj,ix,iy, 0121 & mloc,floc, 0122 & ohat,myThid) 0123 end if 0124 0125 C ==================== reconstruct 5th--order accurate edge values 0126 CALL GAD_PQM_P5E_R(bi,bj,ix,iy, 0127 & mloc,floc, 0128 & edge,myThid) 0129 0130 C ==================== reconstruct coeff. for grid-cell poynomials 0131 CALL GAD_PQM_HAT_R(bi,bj,ix,iy, 0132 & meth, 0133 & mloc,floc, 0134 & edge,ohat, 0135 & fhat,myThid) 0136 0137 C ==================== evaluate integral fluxes on grid-cell edges 0138 CALL GAD_PQM_FLX_R(bi,bj,ix,iy, 0139 & delT,wvel, 0140 & wfac,fhat, 0141 & flux,myThid) 0142 0143 else 0144 0145 do ir = 2, Nr 0146 C ================================== "null" flux on zero transport 0147 flux(ix,iy,ir) = 0. _d 0 0148 end do 0149 0150 end if 0151 0152 end do 0153 end do 0154 0155 return 0156 0157 c end subroutine GAD_PQM_ADV_R 0158 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 |