| 1 |
C $Header: $ |
| 2 |
C $Name: $ |
| 3 |
|
| 4 |
#include "BULK_FORCE_OPTIONS.h" |
| 5 |
|
| 6 |
SUBROUTINE BULKF_SH2RH_AIM(IMODE,NGP,TA,PS,SIG,QA,RH,QSAT,myThid) |
| 7 |
C-- |
| 8 |
C-- SUBROUTINE SHTORH (IMODE,NGP,TA,PS,SIG,QA,RH,QSAT) |
| 9 |
C-- |
| 10 |
C-- Purpose: compute saturation specific humidity and |
| 11 |
C-- relative hum. from specific hum. (or viceversa) |
| 12 |
C-- Input: IMODE : mode of operation |
| 13 |
C-- NGP : no. of grid-points |
| 14 |
C-- TA : abs. temperature |
| 15 |
C-- PS : normalized pressure (= p/1000_hPa) [if SIG < 0] |
| 16 |
C-- : normalized sfc. pres. (= ps/1000_hPa) [if SIG > 0] |
| 17 |
C-- SIG : sigma level |
| 18 |
C-- QA : specific humidity in g/kg [if IMODE = 1] |
| 19 |
C-- RH : relative humidity [if IMODE < 0] |
| 20 |
C-- Output: RH : relative humidity [if IMODE = 1] |
| 21 |
C-- QA : specific humidity in g/kg [if IMODE < 0] |
| 22 |
C-- QSAT : saturation spec. hum. in g/kg |
| 23 |
C-- RH : d.Qsat/d.T in g/kg/K [if IMODE = 2] |
| 24 |
C-- |
| 25 |
|
| 26 |
IMPLICIT NONE |
| 27 |
|
| 28 |
C-- Routine arguments: |
| 29 |
INTEGER IMODE, NGP |
| 30 |
INTEGER myThid |
| 31 |
c _RL TA(NGP), PS(NGP), QA(NGP), RH(NGP), QSAT(NGP) |
| 32 |
_RL TA(NGP), PS(NGP), QSAT(NGP), QA(*), RH(*) |
| 33 |
|
| 34 |
C- jmc: declare all routine arguments: |
| 35 |
_RL SIG |
| 36 |
|
| 37 |
#ifdef ALLOW_BULK_FORCE |
| 38 |
|
| 39 |
C-- Local variables: |
| 40 |
INTEGER J |
| 41 |
|
| 42 |
C- jmc: declare all local variables: |
| 43 |
_RL E0, C1, C2, T0, T1, T2, QS1, QS2 |
| 44 |
_RL sigP, recT, tmpQ |
| 45 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
| 46 |
|
| 47 |
C--- 1. Compute Qsat (g/kg) from T (degK) and normalized pres. P (= p/1000_hPa) |
| 48 |
C If SIG > 0, P = Ps * sigma, otherwise P = Ps(1) = const. |
| 49 |
C |
| 50 |
E0= 6.108 _d -3 |
| 51 |
C1= 17.269 _d 0 |
| 52 |
C2= 21.875 _d 0 |
| 53 |
T0=273.16 _d 0 |
| 54 |
T1= 35.86 _d 0 |
| 55 |
T2= 7.66 _d 0 |
| 56 |
QS1= 622. _d 0 |
| 57 |
QS2= .378 _d 0 |
| 58 |
|
| 59 |
|
| 60 |
IF (IMODE.EQ.2) THEN |
| 61 |
C- Compute Qsat and d.Qsat/d.T : |
| 62 |
DO J=1,NGP |
| 63 |
QSAT(J)=0. |
| 64 |
sigP = PS(1) |
| 65 |
IF (SIG.GT.0.0) sigP=SIG*PS(J) |
| 66 |
IF (TA(J).GE.T0) THEN |
| 67 |
tmpQ = E0*EXP(C1*(TA(J)-T0)/(TA(J)-T1)) |
| 68 |
QSAT(J)= QS1*tmpQ/(sigP-QS2*tmpQ) |
| 69 |
recT = 1. _d 0 / (TA(J)-T1) |
| 70 |
RH(J) = QSAT(J)*C1*(T0-T1)*recT*recT*sigP/(sigP-QS2*tmpQ) |
| 71 |
ELSE IF ( TA(J).GT.T2) THEN |
| 72 |
tmpQ = E0*EXP(C2*(TA(J)-T0)/(TA(J)-T2)) |
| 73 |
QSAT(J)= QS1*tmpQ/(sigP-QS2*tmpQ) |
| 74 |
recT = 1. _d 0 / (TA(J)-T2) |
| 75 |
RH(J) = QSAT(J)*C2*(T0-T2)*recT*recT*sigP/(sigP-QS2*tmpQ) |
| 76 |
ENDIF |
| 77 |
ENDDO |
| 78 |
RETURN |
| 79 |
ENDIF |
| 80 |
|
| 81 |
DO 110 J=1,NGP |
| 82 |
QSAT(J)=0. |
| 83 |
IF (TA(J).GE.T0) THEN |
| 84 |
QSAT(J)=E0*EXP(C1*(TA(J)-T0)/(TA(J)-T1)) |
| 85 |
ELSE IF ( TA(J).GT.T2) THEN |
| 86 |
QSAT(J)=E0*EXP(C2*(TA(J)-T0)/(TA(J)-T2)) |
| 87 |
ENDIF |
| 88 |
110 CONTINUE |
| 89 |
C |
| 90 |
IF (SIG.LE.0.0) THEN |
| 91 |
DO 120 J=1,NGP |
| 92 |
QSAT(J)= QS1*QSAT(J)/( PS(1) - QS2*QSAT(J)) |
| 93 |
120 CONTINUE |
| 94 |
ELSE |
| 95 |
DO 130 J=1,NGP |
| 96 |
QSAT(J)= QS1*QSAT(J)/(SIG*PS(J)- QS2*QSAT(J)) |
| 97 |
130 CONTINUE |
| 98 |
ENDIF |
| 99 |
C |
| 100 |
C--- 2. Compute rel.hum. RH=Q/Qsat (IMODE>0), or Q=RH*Qsat (IMODE<0) |
| 101 |
C |
| 102 |
IF (IMODE.GT.0) THEN |
| 103 |
DO 210 J=1,NGP |
| 104 |
IF(QSAT(J).NE.0.) then |
| 105 |
RH(J)=QA(J)/QSAT(J) |
| 106 |
ELSE |
| 107 |
RH(J)=0. |
| 108 |
ENDIF |
| 109 |
210 CONTINUE |
| 110 |
ELSE IF (IMODE.LT.0) THEN |
| 111 |
DO 220 J=1,NGP |
| 112 |
QA(J)=RH(J)*QSAT(J) |
| 113 |
220 CONTINUE |
| 114 |
ENDIF |
| 115 |
|
| 116 |
#endif /* ALLOW_BULK_FORCE */ |
| 117 |
RETURN |
| 118 |
END |