Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:37:26 UTC

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
b3097ed02d Jean*0001 #include "AIM_OPTIONS.h"
                0002 
                0003 CBOP
                0004 C     !ROUTINE: SUFLUX_LAND
                0005 C     !INTERFACE:
                0006       SUBROUTINE SUFLUX_LAND(
                0007      I                   PSA, FMASK, EMISloc,
                0008      I                   Tsurf, dTskin, SWAV, SSR, SLRD,
e749d70ece Jean*0009      I                   T1, T0, Q0, DENVV,
b3097ed02d Jean*0010      O                   SHF, EVAP, SLRU,
e749d70ece Jean*0011      O                   Shf0, dShf, Evp0, dEvp, Slr0, dSlr, sFlx,
b3097ed02d Jean*0012      O                   TSFC, TSKIN,
                0013      I                   bi,bj,myThid)
                0014 
                0015 C     !DESCRIPTION: \bv
                0016 C     *==========================================================*
                0017 C     | S/R SUFLUX_LAND
                0018 C     | o compute surface flux over land
                0019 C     *==========================================================*
                0020 C     | o contains part of original S/R SUFLUX (Speedy code)
                0021 C     *==========================================================*
                0022 C     \ev
                0023 
                0024 C     !USES:
                0025       IMPLICIT NONE
                0026 
                0027 C     Resolution parameters
                0028 
                0029 C-- size for MITgcm & Physics package :
                0030 #include "AIM_SIZE.h"
                0031 #include "EEPARAMS.h"
                0032 
                0033 C-- Physics package
                0034 #include "AIM_PARAMS.h"
                0035 
                0036 C     Physical constants + functions of sigma and latitude
                0037 #include "com_physcon.h"
                0038 
                0039 C     Surface flux constants
                0040 #include "com_sflcon.h"
                0041 
                0042 C     !INPUT/OUTPUT PARAMETERS:
                0043 C     == Routine Arguments ==
                0044 C--   Input:
                0045 C    PSA    :: norm. surface pressure [p/p0]   (2-dim)
                0046 C    FMASK  :: fractional land-sea mask        (2-dim)
                0047 C    EMISloc:: longwave surface emissivity
                0048 C    Tsurf  :: surface temperature        (2-dim)
                0049 C    dTskin :: temp. correction for daily-cycle heating [K]
                0050 C    SWAV   :: soil wetness availability [0-1] (2-dim)
                0051 C    SSR    :: sfc sw radiation (net flux)     (2-dim)
                0052 C    SLRD   :: sfc lw radiation (downward flux)(2-dim)
e749d70ece Jean*0053 C    T1     :: near-surface air temperature (from Pot.temp)
b3097ed02d Jean*0054 C    T0     :: near-surface air temperature    (2-dim)
                0055 C    Q0     :: near-surface sp. humidity [g/kg](2-dim)
e749d70ece Jean*0056 C    DENVV  :: surface flux (sens,lat.) coeff. (=Rho*|V|) [kg/m2/s]
b3097ed02d Jean*0057 C--   Output:
                0058 C    SHF    :: sensible heat flux              (2-dim)
                0059 C    EVAP   :: evaporation [g/(m^2 s)]         (2-dim)
                0060 C    SLRU   :: sfc lw radiation (upward flux)  (2-dim)
e749d70ece Jean*0061 C    Shf0   :: sensible heat flux over freezing surf.
                0062 C    dShf   :: sensible heat flux derivative relative to surf. temp
b3097ed02d Jean*0063 C    Evp0   :: evaporation computed over freezing surface (Ts=0.oC)
                0064 C    dEvp   :: evaporation derivative relative to surf. temp
                0065 C    Slr0   :: upward long wave radiation over freezing surf.
                0066 C    dSlr   :: upward long wave rad. derivative relative to surf. temp
                0067 C    sFlx   :: net surface flux (+=down) function of surf. temp Ts:
                0068 C              0: Flux(Ts=0.oC) ; 1: Flux(Ts^n) ; 2: d.Flux/d.Ts(Ts^n)
                0069 C    TSFC   :: surface temperature (clim.)     (2-dim)
                0070 C    TSKIN  :: skin surface temperature        (2-dim)
                0071 C--   Input:
                0072 C    bi,bj  :: tile index
                0073 C    myThid :: Thread number for this instance of the routine
                0074 C--
                0075       _RL  PSA(NGP), FMASK(NGP), EMISloc
                0076       _RL  Tsurf(NGP), dTskin(NGP), SWAV(NGP)
                0077       _RL  SSR(NGP), SLRD(NGP) 
e749d70ece Jean*0078       _RL  T1(NGP), T0(NGP), Q0(NGP), DENVV(NGP)
b3097ed02d Jean*0079 
                0080       _RL  SHF(NGP), EVAP(NGP), SLRU(NGP)
e749d70ece Jean*0081       _RL  Shf0(NGP), dShf(NGP), Evp0(NGP), dEvp(NGP)
                0082       _RL  Slr0(NGP), dSlr(NGP), sFlx(NGP,0:2)
b3097ed02d Jean*0083       _RL  TSFC(NGP), TSKIN(NGP)
                0084 
                0085       INTEGER bi,bj,myThid
                0086 CEOP
                0087 
                0088 #ifdef ALLOW_AIM
                0089 
                0090 C-- Local variables:
e749d70ece Jean*0091 C    CDENVV :: surf. heat flux (sens.,lat.) coeff including stability effect
                0092       _RL CDENVV(NGP), RDTH, FSLAND
                0093       _RL Fstb0, dTstb, dFstb
b3097ed02d Jean*0094       _RL QSAT0(NGP,2)
                0095       _RL QDUMMY(1), RDUMMY(1), TS2
                0096       INTEGER J
                0097 
                0098 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0099 
                0100 C     1.5 Define effective skin temperature to compensate for
                0101 C         non-linearity of heat/moisture fluxes during the daily cycle
                0102 
                0103       DO J=1,NGP
                0104         TSKIN(J) = Tsurf(J) + dTskin(J)
                0105         TSFC(J)=273.16 _d 0 + dTskin(J)
                0106       ENDDO
                0107 
                0108 
                0109 C--   2. Computation of fluxes over land and sea
                0110 
e749d70ece Jean*0111 C     2.1 Stability correction
b3097ed02d Jean*0112 
e749d70ece Jean*0113       RDTH = FSTAB/DTHETA
b3097ed02d Jean*0114 
                0115       DO J=1,NGP
e749d70ece Jean*0116         FSLAND=1.+MIN(DTHETA,MAX(-DTHETA,TSKIN(J)-T1(J)))*RDTH
                0117         CDENVV(J)=CHL*DENVV(J)*FSLAND
b3097ed02d Jean*0118       ENDDO
                0119 
e749d70ece Jean*0120       IF ( dTstab.GT.0. _d 0 ) THEN
                0121 C-    account for stability function derivative relative to Tsurf:
                0122 C note: to avoid discontinuity in the derivative (because of min,max), compute
                0123 C   the derivative using the discrete form: F(Ts+dTstab)-F(Ts-dTstab)/2.dTstab
                0124        DO J=1,NGP
                0125         Fstb0 = 1.+MIN(DTHETA,MAX(-DTHETA,TSFC(J) -T1(J)))*RDTH
                0126         Shf0(J) = CHL*DENVV(J)*Fstb0
                0127         dTstb = ( DTHETA+dTstab-ABS(TSKIN(J)-T1(J)) )/dTstab
                0128         dFstb = RDTH*MIN(1. _d 0, MAX(0. _d 0, dTstb*0.5 _d 0))
                0129         dShf(J) = CHL*DENVV(J)*dFstb
                0130        ENDDO
                0131       ENDIF
                0132 
                0133 C     2.2 Evaporation
b3097ed02d Jean*0134 
                0135       CALL SHTORH (2, NGP, TSKIN, PSA, 1. _d 0, QDUMMY, dEvp,
                0136      &                QSAT0(1,1), myThid)
                0137       CALL SHTORH (0, NGP, TSFC, PSA, 1. _d 0, QDUMMY, RDUMMY,
                0138      &                QSAT0(1,2), myThid)
                0139 
e98c426c93 Jean*0140 #ifdef ALLOW_DEW_ON_LAND
                0141 C--   allow negative evaporation (dew):
                0142       IF ( dTstab.GT.0. _d 0 ) THEN
                0143 C-    account for stability function derivative relative to Tsurf:
                0144        DO J=1,NGP
                0145         EVAP(J) = CDENVV(J)*SWAV(J)*(QSAT0(J,1)-Q0(J))
                0146         Evp0(J) =   Shf0(J)*SWAV(J)*(QSAT0(J,2)-Q0(J))
                0147         dEvp(J) = CDENVV(J)*SWAV(J)*dEvp(J)
                0148      &            + dShf(J)*SWAV(J)*(QSAT0(J,1)-Q0(J))
                0149        ENDDO
                0150       ELSE
                0151        DO J=1,NGP
                0152         EVAP(J) = CDENVV(J)*SWAV(J)*(QSAT0(J,1)-Q0(J))
                0153         Evp0(J) = CDENVV(J)*SWAV(J)*(QSAT0(J,2)-Q0(J))
                0154         dEvp(J) = CDENVV(J)*SWAV(J)*dEvp(J)
                0155        ENDDO
                0156       ENDIF
                0157 #else /* ALLOW_DEW_ON_LAND */
                0158 C--   allow only positive evaporation (no dew):
e749d70ece Jean*0159       IF ( dTstab.GT.0. _d 0 ) THEN
                0160 C-    account for stability function derivative relative to Tsurf:
                0161        DO J=1,NGP
                0162         EVAP(J) = CDENVV(J)*SWAV(J)*MAX(0. _d 0,QSAT0(J,1)-Q0(J))
                0163         Evp0(J) =   Shf0(J)*SWAV(J)*MAX(0. _d 0,QSAT0(J,2)-Q0(J))
                0164         dEvp(J) = CDENVV(J)*SWAV(J)*dEvp(J)
                0165      &            + dShf(J)*SWAV(J)*MAX(0. _d 0,QSAT0(J,1)-Q0(J))
                0166        ENDDO
                0167       ELSE
                0168        DO J=1,NGP
b3097ed02d Jean*0169 C       EVAP(J,1) = CDENVV(J,1)*SWAV(J)*MAX(0. _d 0,QSAT0(J,1)-Q0(J))
                0170 c       EVAP(J,1) = CDENVV(J,1)*MAX(0. _d 0,SWAV(J)*QSAT0(J,1)-Q0(J))
                0171 Cjmc: try the other formulation (= described in F.M paper):
                0172         EVAP(J) = CDENVV(J)*SWAV(J)*MAX(0. _d 0,QSAT0(J,1)-Q0(J))
                0173         Evp0(J) = CDENVV(J)*SWAV(J)*MAX(0. _d 0,QSAT0(J,2)-Q0(J))
                0174         dEvp(J) = CDENVV(J)*SWAV(J)*dEvp(J)
e749d70ece Jean*0175        ENDDO
                0176       ENDIF
e98c426c93 Jean*0177 #endif /* ALLOW_DEW_ON_LAND */
e749d70ece Jean*0178 
                0179 C     2.3 Sensible heat flux
                0180 
                0181       IF ( dTstab.GT.0. _d 0 ) THEN
                0182 C-    account for stability function derivative relative to Tsurf:
                0183        DO J=1,NGP
                0184         SHF(J)  = CDENVV(J)*CP*(TSKIN(J)-T0(J))
                0185         Shf0(J) =   Shf0(J)*CP*(TSFC(J) -T0(J))
                0186         dShf(J) = CDENVV(J)*CP
                0187      &            + dShf(J)*CP*(TSKIN(J)-T0(J))
                0188         dShf(J) = MAX( dShf(J), 0. _d 0 )
e98c426c93 Jean*0189 C--   do not allow negative derivative vs Ts of Sensible+Latent H.flux:
                0190 C     a) quiet unrealistic ;
                0191 C     b) garantee positive deriv. of total H.flux (needed for implicit solver)
                0192         dEvp(J) = MAX( dEvp(J), -dShf(J)/ALHC )
e749d70ece Jean*0193        ENDDO
                0194       ELSE
                0195        DO J=1,NGP
                0196         SHF(J)  = CDENVV(J)*CP*(TSKIN(J)-T0(J))
                0197         Shf0(J) = CDENVV(J)*CP*(TSFC(J) -T0(J))
                0198         dShf(J) = CDENVV(J)*CP
                0199        ENDDO
                0200       ENDIF
b3097ed02d Jean*0201 
                0202 C     2.4 Emission of lw radiation from the surface
                0203 
                0204       DO J=1,NGP
                0205         TS2     = TSFC(J)*TSFC(J)
                0206         Slr0(J) = SBC*TS2*TS2
                0207         TS2     = TSKIN(J)*TSKIN(J)
                0208         SLRU(J) = SBC*TS2*TS2
                0209         dSlr(J)  = 4. _d 0 *SBC*TS2*TSKIN(J)
                0210       ENDDO
                0211 
                0212 C--   Compute net surface heat flux (+=down) and its derivative ./. surf. temp.
                0213       DO J=1,NGP
e749d70ece Jean*0214         sFlx(J,0)= ( SSR(J) + SLRD(J) - EMISloc*Slr0(J) )
                0215      &           - ( Shf0(J) + ALHC*Evp0(J) )
                0216         sFlx(J,1)= ( SSR(J) + SLRD(J) - EMISloc*SLRU(J) )
                0217      &           - ( SHF(J)+ ALHC*EVAP(J) )
                0218         sFlx(J,2)=                    - EMISloc*dSlr(J)
                0219      &           - ( dShf(J) + ALHC*dEvp(J) )
b3097ed02d Jean*0220       ENDDO
                0221 
                0222 C--   3. Adjustment of skin temperature and fluxes over land
                0223 C--      based on energy balance (to be implemented)
                0224 C        <= done separately for each surface type (land,sea,sea-ice)
                0225 
                0226 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0227 #endif /* ALLOW_AIM */
                0228 
                0229       RETURN
                0230       END