Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:38:14 UTC

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
88ac963ab8 Jean*0001 #include "BULK_FORCE_OPTIONS.h"
                0002 
                0003       SUBROUTINE BULKF_SH2RH_AIM(IMODE,NGP,TA,PS,SIG,QA,RH,QSAT,myThid)
                0004 C--
                0005 C--   SUBROUTINE SHTORH (IMODE,NGP,TA,PS,SIG,QA,RH,QSAT)
                0006 C--
                0007 C--   Purpose: compute saturation specific humidity and
                0008 C--            relative hum. from specific hum. (or viceversa)
                0009 C--   Input:   IMODE  : mode of operation
                0010 C--            NGP    : no. of grid-points
                0011 C--            TA     : abs. temperature
                0012 C--            PS     : normalized pressure   (=  p/1000_hPa) [if SIG < 0]
                0013 C--                   : normalized sfc. pres. (= ps/1000_hPa) [if SIG > 0]
                0014 C--            SIG    : sigma level
                0015 C--            QA     : specific humidity in g/kg [if IMODE = 1]
                0016 C--            RH     : relative humidity         [if IMODE < 0]
                0017 C--   Output:  RH     : relative humidity         [if IMODE = 1]
                0018 C--            QA     : specific humidity in g/kg [if IMODE < 0]
                0019 C--            QSAT   : saturation spec. hum. in g/kg
                0020 C--            RH     : d.Qsat/d.T  in g/kg/K     [if IMODE = 2]
                0021 C--
                0022 
                0023       IMPLICIT NONE
                0024 
                0025 C-- Routine arguments:
                0026       INTEGER IMODE, NGP
                0027       INTEGER  myThid
                0028 c     _RL TA(NGP), PS(NGP), QA(NGP), RH(NGP), QSAT(NGP)
                0029       _RL TA(NGP), PS(NGP), QSAT(NGP), QA(*), RH(*)
                0030 
                0031 C- jmc: declare all routine arguments:
                0032       _RL SIG
                0033 
                0034 #ifdef ALLOW_BULK_FORCE
                0035 
                0036 C-- Local variables:
                0037       INTEGER  J
                0038 
                0039 C- jmc: declare all local variables:
                0040       _RL E0, C1, C2, T0, T1, T2, QS1, QS2
                0041       _RL sigP, recT, tmpQ
                0042 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0043 
                0044 C---  1. Compute Qsat (g/kg) from T (degK) and normalized pres. P (= p/1000_hPa)
                0045 C        If SIG > 0, P = Ps * sigma, otherwise P = Ps(1) = const.
                0046 C
                0047       E0=  6.108 _d -3
                0048       C1= 17.269 _d 0
                0049       C2= 21.875 _d 0
                0050       T0=273.16 _d 0
                0051       T1= 35.86 _d 0
                0052       T2=  7.66 _d 0
                0053       QS1= 622. _d 0
                0054       QS2= .378 _d 0
                0055 
                0056 
                0057       IF (IMODE.EQ.2) THEN
                0058 C-    Compute Qsat and d.Qsat/d.T :
                0059         DO J=1,NGP
                0060          QSAT(J)=0.
                0061          sigP = PS(1)
                0062          IF (SIG.GT.0.0) sigP=SIG*PS(J)
                0063          IF (TA(J).GE.T0) THEN
                0064           tmpQ   = E0*EXP(C1*(TA(J)-T0)/(TA(J)-T1))
                0065           QSAT(J)= QS1*tmpQ/(sigP-QS2*tmpQ)
                0066           recT   = 1. _d 0 / (TA(J)-T1)
                0067           RH(J)  = QSAT(J)*C1*(T0-T1)*recT*recT*sigP/(sigP-QS2*tmpQ)
                0068          ELSE IF ( TA(J).GT.T2) THEN
                0069           tmpQ   = E0*EXP(C2*(TA(J)-T0)/(TA(J)-T2))
                0070           QSAT(J)= QS1*tmpQ/(sigP-QS2*tmpQ)
                0071           recT   = 1. _d 0 / (TA(J)-T2)
                0072           RH(J)  = QSAT(J)*C2*(T0-T2)*recT*recT*sigP/(sigP-QS2*tmpQ)
                0073          ENDIF
                0074         ENDDO
                0075         RETURN
                0076       ENDIF
                0077 
                0078       DO 110 J=1,NGP
                0079         QSAT(J)=0.
                0080         IF (TA(J).GE.T0) THEN
                0081           QSAT(J)=E0*EXP(C1*(TA(J)-T0)/(TA(J)-T1))
                0082         ELSE IF ( TA(J).GT.T2) THEN
                0083           QSAT(J)=E0*EXP(C2*(TA(J)-T0)/(TA(J)-T2))
                0084         ENDIF
                0085   110 CONTINUE
                0086 C
                0087       IF (SIG.LE.0.0) THEN
                0088         DO 120 J=1,NGP
                0089           QSAT(J)= QS1*QSAT(J)/( PS(1)   - QS2*QSAT(J))
                0090   120   CONTINUE
                0091       ELSE
                0092         DO 130 J=1,NGP
                0093           QSAT(J)= QS1*QSAT(J)/(SIG*PS(J)- QS2*QSAT(J))
                0094   130   CONTINUE
                0095       ENDIF
                0096 C
                0097 C---  2. Compute rel.hum. RH=Q/Qsat (IMODE>0), or Q=RH*Qsat (IMODE<0)
                0098 C
                0099       IF (IMODE.GT.0) THEN
                0100         DO 210 J=1,NGP
                0101           IF(QSAT(J).NE.0.) then
                0102             RH(J)=QA(J)/QSAT(J)
                0103           ELSE
                0104             RH(J)=0.
                0105           ENDIF
                0106   210   CONTINUE
                0107       ELSE IF (IMODE.LT.0) THEN
                0108         DO 220 J=1,NGP
                0109           QA(J)=RH(J)*QSAT(J)
                0110   220   CONTINUE
                0111       ENDIF
                0112 
                0113 #endif /* ALLOW_BULK_FORCE */
                0114       RETURN
                0115       END