!-------------------------------------- LICENCE BEGIN ------------------------------------
!Environment Canada - Atmospheric Science and Technology License/Disclaimer,
! version 3; Last Modified: May 7, 2008.
!This is free but copyrighted software; you can use/redistribute/modify it under the terms
!of the Environment Canada - Atmospheric Science and Technology License/Disclaimer
!version 3 or (at your option) any later version that should be found at:
!http://collaboration.cmc.ec.gc.ca/science/rpn.comm/license.html
!
!This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
!without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
!See the above mentioned License/Disclaimer for more details.
!You should have received a copy of the License/Disclaimer along with this software;
!if not, you can write to: EC-RPN COMM Group, 2121 TransCanada, suite 500, Dorval (Quebec),
!CANADA, H9P 1J3; or send e-mail to service.rpn@ec.gc.ca
!-------------------------------------- LICENCE END --------------------------------------
***S/P MCICA_CLD_GENERATOR - PART OF THE ISCCP CLOUD SIMULATOR PACKAGE
*
SUBROUTINE MCICA_CLD_GENERATOR(ZF, AVG_CF, AVG_QLW,AVG_QIW, 1
1 RLC_CF, RLC_CW,SIGMA_QCW,
2 ILG, IL1, IL2, LAY,
3 ICEWCIN_S,LIQWCIN_S, NCLDY)
#include "impnone.cdk"
#include "mcica.cdk"
#include "consphy.cdk"
!#include "xcwdata.cdk"
! INPUT DATA
!
REAL, INTENT(IN) ::
1 ZF(ILG,LAY), ! Full-level (layer midpoint) altitude (km)
2 AVG_CF(ILG,LAY), ! Cloud fraction for each layer
3 AVG_QLW(ILG,LAY), ! Cloud mean liquid water mixing ratio (g/m^3)
4 AVG_QIW(ILG,LAY), ! Cloud mean ice mixing ratio (g/m^3)
5 SIGMA_QCW(ILG,LAY), ! Normalized cloud condensate std. dev.
6 RLC_CF(ILG,LAY), ! Cloud fraction decorrelation length
7 RLC_CW(ILG,LAY) ! Cloud condensate decorrelation length
INTEGER, INTENT(IN) ::
1 ILG,
2 IL1,
3 IL2,
4 LAY
!
! OUTPUT DATA
!
REAL,INTENT(OUT) ::
1 ICEWCIN_S(ILG,LAY,NX_LOC), ! Column ice water mixing ratio profile (g/m^3)
2 LIQWCIN_S(ILG,LAY,NX_LOC) ! Column liquid water mixing ratio profile (g/m^3)
INTEGER, INTENT(OUT) ::
1 NCLDY(ILG) ! Number of cloudy subcolumns
*Author
* Jason Cole, MSC/Cloud Physics (Summer 2005)
*Revisions
* 001 ...
*Object
* Generate a NX column by LAY layer cloud field. Input profiles
* must be from top of the model to the bottom as is the output.
* A Method and model evaluation is described in:
* Raisanen et al, 2004, Stochastic generation of subgrid-scale
* cloudy columns for large-scale models, Quart. J. Roy. Meteorol.
* Soc., 2004 , 130 , 2047-2068
*Arguments
* - Input -
* ZF (ILG,LAY) Full-level (layer midpoint) altitude (km)
* AVG_CF (ILG,LAY) Cloud fraction for each layer
* AVG_QLW (ILG,LAY) Cloud mean liquid water mixing ratio (g/m^3)
* AVG_QIW (ILG,LAY) Cloud mean ice mixing ratio (g/m^3)
* SIGMA_QCW(ILG,LAY) Normalized cloud condensate std. dev.
* RLC_CF (ILG,LAY) Cloud fraction decorrelation length
* RLC_CW (ILG,LAY) Cloud condensate decorrelation length
* ILG, IL1, IL2 Horizontal dimension
* LAY Number of layers in GCM
* NX_LOC Number of columns to generate
* - Ouput -
* ICEWCIN_S(ILG,LAY,NX_LOC) Column ice water mixing ratio profile (g/m^3)
* LIQWCIN_S(ILG,LAY,NX_LOC) Column liquid water mixing ratio profile (g/m^3)
* NCLDY (ILG) Number of cloudy subcolumns
*Implicites
REAL, PARAMETER ::
1 CUT = 0.001
!
! LOCAL DATA
!
INTEGER ::
1 IL, ! Counter over ILG GCM columns
2 I, ! Counter
2 K, ! Counter over LAY vertical layers
3 K_TOP, ! Index of top most cloud layer
4 K_BASE, ! Index of lowest cloud layer
5 IND1, ! Index in variability calculation
6 IND2 ! Index in variability calculation
REAL ::
1 ALPHA(ILG,LAY), ! Fraction of maximum/random cloud overlap
2 RCORR(ILG,LAY), ! Fraction of maximum/random cloud condensate overlap
5 RIND1, ! Real index in variability calculation
6 RIND2, ! Real index in variability calculation
7 ZCW ! Ratio of cloud condensate miximg ratio for this cell to its layer cloud-mean value
REAL :: ! Random number vectors
1 X(ILG),
2 Y(ILG),
3 X1(ILG),
4 Y1(ILG),
5 X2(ILG),
6 Y2(ILG)
INTEGER ::
1 I_LOC(ILG) ! Place the new subcolumns into the arrays starting from the front
LOGICAL ::
1 L_TOP,
2 L_BOT,
3 L_CLD(ILG),
4 L_CLDSUB(ILG)
! Initialize the arrays
DO IL = IL1, IL2
I_LOC(IL) = 1
L_CLD(IL) = .FALSE.
END DO
! Find uppermost cloudy layer
L_TOP = .FALSE.
K_TOP = 0
DO K=1,LAY
DO IL = IL1, IL2
K_TOP = K
IF (AVG_CF(IL,K) .GT. CUT) L_TOP = .TRUE.
END DO ! IL
IF (L_TOP) EXIT
END DO ! K
! If no cloudy layers in any GCM column in this group, exit
IF (K_TOP .EQ. 0) THEN
RETURN
END IF
! Find lowermost cloudy layer
L_BOT = .FALSE.
K_BASE = 0
DO K=LAY,1,-1
DO IL = IL1, IL2
K_BASE = K
IF (AVG_CF(IL,K) .GT. CUT) L_BOT = .TRUE.
END DO ! IL
IF (L_BOT) EXIT
END DO ! K
! Calculate overlap factors ALPHA for cloud fraction and RCORR for cloud
! condensate based on layer midpoint distances and decorrelation depths
DO K=1, LAY !K_TOP,K_BASE-1
DO IL = IL1, IL2
! IF (RLC_CF(IL,K) .GT. 0.0) THEN
! ALPHA(IL,K) = EXP(-(ZF(IL,K) - ZF(IL,K+1)) / RLC_CF(IL,K))
ALPHA(IL,K) = -(ZF(IL,K) - ZF(IL,K+1)) / RLC_CF(IL,K)
! ELSE
! ALPHA(IL,K) = 0.0
! END IF
! IF (RLC_CW(IL,K) .GT. 0.0) THEN
! RCORR(IL,K) = EXP(-(ZF(IL,K) - ZF(IL,K+1)) / RLC_CW(IL,K))
RCORR(IL,K) = -(ZF(IL,K) - ZF(IL,K+1)) / RLC_CW(IL,K)
! ELSE
! RCORR(IL,K) = 0.0
! END IF
END DO ! IL
END DO ! K
CALL VSEXP(ALPHA,ALPHA,LAY*(IL2-IL1+1))
CALL VSEXP(RCORR,RCORR,LAY*(IL2-IL1+1))
DO I=1,NX_LOC
! Generate all subcolumns for latitude chain
DO K = K_TOP, K_BASE
IF (K .EQ. K_TOP) THEN
CALL RANDOM_NUMBER(X)
CALL RANDOM_NUMBER(Y)
END IF
CALL RANDOM_NUMBER(X1)
CALL RANDOM_NUMBER(Y1)
CALL RANDOM_NUMBER(X2)
CALL RANDOM_NUMBER(Y2)
DO IL = IL1, IL2
! Maximum-random overlap
IF (LMAXRAN) THEN
IF (X(IL) .LT. 1.0-AVG_CF(IL,K-1)) THEN !It is clear above
X(IL) = X1(IL) * (1.0 - AVG_CF(IL,K-1))
Y(IL) = Y1(IL)
END IF
! Generalized overlap
ELSE
IF (X1(IL) .GT. ALPHA(IL,K-1)) X(IL) = X2(IL)
IF (Y1(IL) .GT. RCORR(IL,K-1)) Y(IL) = Y2(IL)
END IF
! Treatment of cloudy cells
IF (X(IL) .GT. 1.0-AVG_CF(IL,K)) THEN ! Generate cloud in this layer
IF (LPPH) THEN ! Homogeneous clouds
ZCW = 1.0
ELSE
! Horizontally variable clouds:
! Determine ZCW = ratio of cloud condensate miximg ratio QC for this cell to
! its mean value for all cloudy cells in this layer.
! Use bilinear interpolation of ZCW tabulated in array XCW as a function of
! * cumulative probability Y
! * relative standard deviation SIGMA
! Take care that the definition of RIND2 is consistent with subroutine
! TABULATE_XCW
RIND1 = Y(IL) * (N1 - 1) + 1.0
IND1 = MAX(1, MIN(INT(RIND1), N1-1))
RIND1 = RIND1 - IND1
RIND2 = 40.0 * SIGMA_QCW(IL,K) - 3.0
IND2 = MAX(1, MIN(INT(RIND2), N2-1))
RIND2 = RIND2 - IND2
ZCW = (1.0-RIND1) * (1.0-RIND2) * XCW(IND1,IND2)
1 + (1.0-RIND1) * RIND2 * XCW(IND1,IND2+1)
2 + RIND1 * (1.0-RIND2) * XCW(IND1+1,IND2)
3 + RIND1 * RIND2 * XCW(IND1+1,IND2+1)
END IF
! A horizontally constant IWC/LWC ratio is assumed for each layer so far
L_CLD(IL) = .TRUE.
LIQWCIN_S(IL,K,I_LOC(IL)) = ZCW * AVG_QLW(IL,K)
ICEWCIN_S(IL,K,I_LOC(IL)) = ZCW * AVG_QIW(IL,K)
END IF
END DO ! IL
END DO ! K
! Need to check if a cloudy subcolumn was generated
DO IL = IL1, IL2
IF (L_CLD(IL)) THEN
I_LOC(IL) = I_LOC(IL) + 1
L_CLD(IL) = .FALSE.
END IF
END DO
END DO ! I
! Record the number of cloudy subcolumns generated
DO IL = IL1, IL2
NCLDY(IL) = I_LOC(IL) - 1
END DO ! IL
RETURN
END SUBROUTINE MCICA_CLD_GENERATOR
SUBROUTINE PREP_MCICA(RLC_CF, RLC_CW, SIGMA_QCW, NUAGE, ILG, IL1, 1
1 IL2, LAY)
! --------------------------------------------------------------------
! Subroutine to define the horizontal variability and overlap for the
! stochastic cloud generator.
! --------------------------------------------------------------------
IMPLICIT NONE
!
! INPUT DATA
!
INTEGER, INTENT(IN) ::
1 ILG,
2 IL1,
3 IL2,
4 LAY
REAL, INTENT(IN) ::
1 NUAGE(ILG,LAY) ! Cloud amount
!
! OUTPUT DATA
!
REAL, INTENT(OUT) ::
1 RLC_CF(ILG,LAY), ! Cloud fraction decorrelation length (km)
2 RLC_CW(ILG,LAY), ! Cloud condensate decorrelation length (km)
3 SIGMA_QCW(ILG,LAY) ! Normalized standard deviation of cloud condensate (unitless)
!
! LOCAL DATA
!
INTEGER ::
1 IL,
2 KK
REAL ::
1 ANU
! ZERO OUT THE OUTPUT ARRAYS
DO KK = 1, LAY
DO IL = IL1, IL2
SIGMA_QCW(IL,KK) = 0.0
RLC_CF(IL,KK) = 0.0
RLC_CW(IL,KK) = 0.0
END DO ! IL
END DO ! KK
! FOR NOW SET THE DECORRELATION LENGTHS TO BE REASONABLE ESTIMATES
! UPDATE WHEN HAVE SOME CLUE ABOUT HOW TO DO THIS BASED ON LARGE-SCALE
! VARIABLES
DO KK = 1, LAY
DO IL = IL1, IL2
RLC_CF(IL,KK) = 2.0
RLC_CW(IL,KK) = 1.0
END DO
END DO
! DEFINE THE HORIZONTAL VARIABILITY USING METHOD CURRENTLY USED IN GEM
! WILL CHANGE SOON
DO KK = 1, LAY
DO IL = IL1, IL2
IF (NUAGE(IL,KK) .LE. 0.9) THEN
ANU = 1.0
ELSEIF (NUAGE(IL,KK) .GT. 0.9 .AND.
1 NUAGE(IL,KK) .LT. 1.0) THEN
ANU = 2.0
ELSEIF (NUAGE(IL,KK) .GT. 0.0) THEN
ANU = 4.0
ENDIF
! THE NORMALIZED VARIABILITY (STD DEV/MEAN) CAN BE APPROXIMATED AS
! 1.0/SQRT(ANU)
IF (NUAGE(IL,KK) .GT. 0.0) THEN
SIGMA_QCW(IL,KK) = 1.0/SQRT(ANU)
ELSE
SIGMA_QCW(IL,KK) = 0.0
END IF
END DO ! IL
END DO ! KK
RETURN
END SUBROUTINE PREP_MCICA
SUBROUTINE McICA_CLD_GEN(NUAGE, LIQWCIN, ICEWCIN, RLC_CF, RLC_CW, 1,1
1 SIGMA_QCW, T, S, PS, ILG, IL1, IL2, LAY,
2 NCLDY, LIQWCIN_S, ICEWCIN_S, CLDTOT)
! --------------------------------------------------------------------
! Driver to call stochastic cloud generator and produce diagnostic cloud
! properties that are consistent with McICA radiative transfer routine.
! --------------------------------------------------------------------
IMPLICIT NONE
! EXTERNAL MCICA_CLD_GENERATOR
#include "mcica.cdk"
#include "consphy.cdk"
! Note: LAY => Number of layers
! Note: NX_LOC => Number of subcolumns to generate
!
! PARAMETER
!
REAL, PARAMETER ::
1 M2KM = 1.0/1000.0, ! Convert meters to kilometers
2 CUT = 0.001
!
! INPUT DATA
!
REAL, INTENT(IN) ::
1 NUAGE(ILG,LAY), ! Column cloud fraction
2 ICEWCIN(ILG,LAY), ! Column in-cloud ice water mixing ratio profile (kg/kg)
3 LIQWCIN(ILG,LAY) ! Column in-cloud liquid water mixing ratio profile (kg/kg)
INTEGER, INTENT(IN) :: ! Counters and array sizes
1 ILG,
2 IL1,
3 IL2,
4 LAY
REAL, INTENT(IN) ::
1 RLC_CF(ILG,LAY), ! Cloud fraction decorrelation length (km)
2 RLC_CW(ILG,LAY), ! Cloud condensate decorrelation length (km)
3 SIGMA_QCW(ILG,LAY), ! Normalized standard deviation of cloud condensate (unitless)
4 T(ILG,LAY), ! Column temperature at layer midpoint (K)
5 PS(ILG), ! Surface pressure (Pa)
6 S(ILG,LAY) ! Column hybrid coord at layer midpoint (unitless)
!
! OUTPUT DATA
!
REAL, INTENT(OUT) ::
1 CLDTOT(ILG), ! McICA vertical projected total cloud fraction (unitless)
2 ICEWCIN_S(ILG,LAY,NX_LOC), ! Column ice water mixing ratio profile (g/m^3)
3 LIQWCIN_S(ILG,LAY,NX_LOC) ! Column liquid water mixing ratio profile (g/m^3)
INTEGER,INTENT(OUT) ::
1 NCLDY(ILG) ! Number of cloudy subcolumns
!
! LOCAL DATA
!
INTEGER ::
1 IL, ! Counter over GCM columns
2 II, ! Counter over NX subcolumns
3 KK ! Counter over lay vertical layers
REAL ::
1 RHO ! Density of air (g/m^3)
REAL ::
1 P, ! GCM column pressure at layer midpoint (Pa)
2 ROG ! Total gas constant/gravity
REAL ::
1 DMULT,
2 Q_TOT
REAL ::
1 X1(ILG)
REAL ::
1 QI_PROF(ILG,LAY),
2 QC_PROF(ILG,LAY),
3 ZM(ILG,LAY),
4 DZ(ILG,LAY)
INTEGER ::
1 K,
2 IND1
! Zero out fields
DO IL = IL1,IL2
CLDTOT(IL) = 0.0
NCLDY(IL) = 0
END DO ! IL
DO KK = 1, LAY
DO IL = IL1,IL2
QC_PROF(IL,KK) = 0.0
QI_PROF(IL,KK) = 0.0
END DO ! IL
END DO ! KK
DO II = 1, NX_LOC
DO KK = 1 , LAY
DO IL = IL1,IL2
ICEWCIN_S(IL,KK,II) = 0.0
LIQWCIN_S(IL,KK,II) = 0.0
END DO
END DO
END DO
! Compute the heights of mid-layers
ROG=RGASD/GRAV
DO KK = 1, LAY
DO IL = IL1, IL2
P = S(IL,KK)*PS(IL)
ZM(IL,KK) = ROG*T(IL,KK)*LOG(PS(IL)/P)*M2KM
END DO
END DO
! Convert the cloud condensate from kg/kg to g/m^3 (is the cloud condensate cloud or domain mean)?
DO KK = 1, LAY
DO IL = IL1, IL2
P = S(IL,KK)*PS(IL)
! Compute layer height
IF (NUAGE(IL,KK) .GT. CUT) THEN
! RHO in kg/m^3
RHO = 1000.0*P/(RGASD*T(IL,KK))
DMULT = RHO
QI_PROF(IL,KK) = ICEWCIN(IL,KK)*DMULT
QC_PROF(IL,KK) = LIQWCIN(IL,KK)*DMULT
END IF
END DO ! IL
END DO ! KK
! Call cloud generator
CALL MCICA_CLD_GENERATOR
(ZM, NUAGE, QC_PROF, QI_PROF,
1 RLC_CF, RLC_CW, SIGMA_QCW,
2 ILG, IL1, IL2, LAY,
4 ICEWCIN_S,LIQWCIN_S, NCLDY)
DO IL = 1, ILG
CLDTOT(IL) = REAL(NCLDY(IL))/REAL(NX_LOC)
END DO !IL
RETURN
END SUBROUTINE McICA_CLD_GEN