source: LMDZ6/trunk/libf/phylmd/lmdz_cloud_optics_prop.F90 @ 4707

Last change on this file since 4707 was 4707, checked in by Laurent Fairhead, 12 months ago

Continuing on from the morning's poihl workshop: getting rid of the includes in the routine

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 27.1 KB
Line 
1! $Id: lmdz_cloud_optics_prop.F90 4707 2023-09-25 18:19:15Z fairhead $
2MODULE lmdz_cloud_optics_prop
3
4CONTAINS
5
6SUBROUTINE cloud_optics_prop(klon, klev, paprs, pplay, temp, radocond, picefra, pclc, &
7    pcltau, pclemi, pch, pcl, pcm, pct, radocondwp, xflwp, xfiwp, xflwc, xfiwc, &
8    mass_solu_aero, mass_solu_aero_pi, pcldtaupi, distcltop, temp_cltop, re, fl, reliq, reice, &
9    reliq_pi, reice_pi, scdnc, cldncl, reffclwtop, lcc, reffclws, &
10    reffclwc, cldnvi, lcc3d, lcc3dcon, lcc3dstra, icc3dcon, icc3dstra,  &
11    icefrac_optics, dNovrN, ptconv,rnebcon, ccwcon)
12
13  USE lmdz_cloud_optics_prop_ini , ONLY : flag_aerosol, ok_cdnc
14  USE lmdz_cloud_optics_prop_ini , ONLY : lunout
15  USE lmdz_cloud_optics_prop_ini , ONLY : bl95_b0, bl95_b1
16  USE lmdz_cloud_optics_prop_ini , ONLY : latitude_deg
17  USE lmdz_cloud_optics_prop_ini , ONLY : iflag_t_glace
18  USE lmdz_cloud_optics_prop_ini , ONLY : cdnc_max, cdnc_max_m3
19  USE lmdz_cloud_optics_prop_ini , ONLY : cdnc_min, cdnc_min_m3
20  USE lmdz_cloud_optics_prop_ini , ONLY : thres_tau, thres_neb
21  USE lmdz_cloud_optics_prop_ini , ONLY : prmhc, prlmc
22  USE lmdz_cloud_optics_prop_ini , ONLY : coef_froi, coef_chau
23  USE lmdz_cloud_optics_prop_ini , ONLY : seuil_neb
24  USE lmdz_cloud_optics_prop_ini , ONLY : t_glace_min_old, t_glace_max_old
25  USE lmdz_cloud_optics_prop_ini , ONLY : k_ice0, df
26  USE lmdz_cloud_optics_prop_ini , ONLY : rg, rd, rpi
27  USE lmdz_cloud_optics_prop_ini , ONLY : rad_chau1, rad_chau2, rad_froid, iflag_rei
28  USE lmdz_cloud_optics_prop_ini , ONLY : ok_icefra_lscp, rei_max, rei_min
29  USE lmdz_cloud_optics_prop_ini , ONLY : zepsec, novlp, iflag_ice_thermo, ok_new_lscp
30 
31  USE icefrac_lsc_mod ! computes ice fraction (JBM 3/14)
32  USE lmdz_lscp_tools, only: icefrac_lscp
33
34
35
36  IMPLICIT NONE
37  ! ======================================================================
38  ! Authors: Z.X. Li (LMD/CNRS) date: 19930910
39  !          O.Boucher (LMD/CNRS) mise a jour en 201212
40  !          I. Musat (LMD/CNRS) : prise en compte de la meme hypothese
41  !                              de recouvrement pour les nuages que pour
42  !                              le rayonnement rrtm via le parametre
43  !                               novlp de radopt.h : 20160721
44  !          L.Fairheard, E.Vignon, JB Madeleine, L. Raillard, A. Idelkadi
45  !          M. Coulon-Decorzens: replayisation of the routine + cleaning
46  !                               and commentaries
47  !
48  ! Aim: compute condensate optical properties,
49  !      cloud optical depth and emissivity
50  ! ======================================================================
51 
52  ! List of arguments
53  !------------------
54
55  ! input:
56  INTEGER, INTENT(IN) :: klon, klev      ! number of horizontal and vertical grid points
57  REAL, INTENT(IN) :: paprs(klon, klev+1)! pressure at bottom interfaces [Pa]
58  REAL, INTENT(IN) :: pplay(klon, klev)  ! pressure at the middle of layers [Pa]
59  REAL, INTENT(IN) :: temp(klon, klev)   ! temperature [K]
60  REAL, INTENT(IN) :: radocond(klon, klev) ! cloud condensed water seen by radiation [kg/kg]
61  REAL, INTENT(IN) :: picefra(klon,klev) ! ice fraction in clouds from large scale condensation scheme [-]
62  REAL, INTENT(IN) :: rnebcon(klon,klev) ! convection cloud fraction [-]
63  REAL, INTENT(IN) :: ccwcon(klon,klev)  ! condensed water from deep convection [kg/kg]
64  ! jq for the aerosol indirect effect
65  ! jq introduced by Johannes Quaas (quaas@lmd.jussieu.fr), 27/11/2003
66  REAL, INTENT(IN) :: mass_solu_aero(klon, klev)    ! total mass concentration for all soluble aerosols [ug m-3]
67  REAL, INTENT(IN) :: mass_solu_aero_pi(klon, klev) ! - (pre-industrial value)
68  REAL, INTENT(IN)  :: dNovrN(klon)         ! enhancement factor for cdnc
69  REAL, INTENT(OUT) :: distcltop(klon,klev) ! distance from large scale cloud top [m]
70  REAL, INTENT(OUT) :: temp_cltop(klon,klev)!temperature at large scale cloud top [K]
71
72  LOGICAL, INTENT(IN) :: ptconv(klon, klev) ! flag for grid points affected by deep convection
73
74  ! inout:
75  REAL, INTENT(INOUT) :: pclc(klon, klev) ! cloud fraction for radiation [-]
76
77  ! out:
78  REAL, INTENT(OUT) :: pct(klon)      ! 2D total cloud cover [-]
79  REAL, INTENT(OUT) :: pcl(klon)      ! 2D low cloud cover [-]
80  REAL, INTENT(OUT) :: pcm(klon)      ! 2D mid cloud cover [-]
81  REAL, INTENT(OUT) :: pch(klon)      ! 2D high cloud cover [-]
82  REAL, INTENT(OUT) :: radocondwp(klon) ! total condensed water path (seen by radiation) [kg/m2]
83  REAL, INTENT(OUT) :: xflwp(klon)    ! liquid water path (seen by radiation) [kg/m2]
84  REAL, INTENT(OUT) :: xfiwp(klon)    ! ice water path (seen by radiation) [kg/m2]
85  REAL, INTENT(OUT) :: xflwc(klon, klev) ! liquid water content seen by radiation [kg/kg]
86  REAL, INTENT(OUT) :: xfiwc(klon, klev) ! ice water content seen by radiation [kg/kg]
87  REAL, INTENT(OUT) :: re(klon, klev) ! cloud droplet effective radius multiplied by fl
88  REAL, INTENT(OUT) :: fl(klon, klev) ! xliq * rneb, denominator to re; fraction of liquid water clouds
89                                      ! introduced to avoid problems in the averaging of the output
90                                      ! water clouds within a grid cell
91
92  REAL, INTENT(OUT) :: pcltau(klon, klev) ! cloud optical depth [m]
93  REAL, INTENT(OUT) :: pclemi(klon, klev) ! cloud emissivity [-]
94  REAL, INTENT(OUT) :: pcldtaupi(klon, klev) ! pre-industrial value of cloud optical thickness, ie.
95                                             ! values of optical thickness that does not account
96                                             ! for aerosol effects on cloud droplet radius [m]
97
98  REAL, INTENT(OUT) :: reliq(klon, klev)   ! liquid droplet effective radius [m]
99  REAL, INTENT(OUT) :: reice(klon, klev)   ! ice effective radius [m]
100  REAL, INTENT(OUT) :: reliq_pi(klon, klev)! liquid droplet effective radius [m], pre-industrial
101  REAL, INTENT(OUT) :: reice_pi(klon, klev)! ice effective radius [m], pre-industrial
102  REAL, INTENT(OUT) :: scdnc(klon, klev)   ! cloud droplet number concentration, mean over the whole mesh [m-3]
103  REAL, INTENT(OUT) :: cldncl(klon)        ! cloud droplet number concentration at top of cloud [m-3]
104  REAL, INTENT(OUT) :: reffclwtop(klon)    ! effective radius of cloud droplet at top of cloud [m]
105  REAL, INTENT(OUT) :: lcc(klon)           ! liquid Cloud Content at top of cloud [kg/kg]
106  REAL, INTENT(OUT) :: reffclws(klon, klev)! stratiform cloud droplet effective radius
107  REAL, INTENT(OUT) :: reffclwc(klon, klev)! convective cloud droplet effective radius
108  REAL, INTENT(OUT) :: cldnvi(klon)        ! column Integrated cloud droplet Number [/m2]
109  REAL, INTENT(OUT) :: lcc3d(klon, klev)   ! cloud fraction for liquid part only [-]
110  REAL, INTENT(OUT) :: lcc3dcon(klon, klev)! cloud fraction for liquid part only, convective clouds [-]
111  REAL, INTENT(OUT) :: lcc3dstra(klon, klev)!cloud fraction for liquid part only, stratiform clouds [-]
112  REAL, INTENT(OUT) :: icc3dcon(klon, klev)! cloud fraction for liquid part only, convective clouds [-]
113  REAL, INTENT(OUT) :: icc3dstra(klon, klev)! cloud fraction for ice part only, stratiform clouds [-]
114  REAL, INTENT(OUT) :: icefrac_optics(klon, klev)! ice fraction in clouds seen by radiation [-]
115
116  ! Local variables
117  !----------------
118
119  LOGICAL, SAVE :: first = .TRUE.
120  !$OMP THREADPRIVATE(FIRST)
121  INTEGER flag_max
122
123  ! threshold PARAMETERs
124  REAL phase3d(klon, klev)
125  REAL tcc(klon), ftmp(klon), lcc_integrat(klon), height(klon)
126  LOGICAL lo
127  INTEGER i, k
128  REAL radius
129
130
131  REAL rel, tc, rei, iwc, dei, deimin, deimax
132  REAL k_ice
133
134  ! jq for the aerosol indirect effect
135  ! jq introduced by Johannes Quaas (quaas@lmd.jussieu.fr), 27/11/2003
136  REAL cdnc(klon, klev) ! cloud droplet number concentration [m-3]
137  REAL cdnc_pi(klon, klev) ! cloud droplet number concentration [m-3] (pi value)
138  REAL re_pi(klon, klev) ! cloud droplet effective radius [um] (pi value)
139
140  ! IM cf. CR:parametres supplementaires
141  REAL dzfice(klon,klev)
142  REAL zclear(klon)
143  REAL zcloud(klon)
144  REAL zcloudh(klon)
145  REAL zcloudm(klon)
146  REAL zcloudl(klon)
147  REAL rhodz(klon, klev) !--rho*dz pour la couche
148  REAL zrho(klon, klev) !--rho pour la couche
149  REAL dh(klon, klev) !--dz pour la couche
150  REAL rad_chaud(klon, klev) !--rayon pour les nuages chauds
151  REAL rad_chaud_pi(klon, klev) !--rayon pour les nuages chauds pre-industriels
152  REAL zflwp_var, zfiwp_var
153  REAL d_rei_dt
154
155
156  ! FH : 2011/05/24
157  ! rei = ( rei_max - rei_min ) * T(°C) / 81.4 + rei_max
158  ! to be used for a temperature in celcius T(°C) < 0
159  ! rei=rei_min for T(°C) < -81.4
160  ! Calcul de la pente de la relation entre rayon effective des cristaux
161  ! et la température Pour retrouver les résultats numériques de la version d'origine,
162  ! on impose 0.71 quand on est proche de 0.71
163  d_rei_dt = (rei_max-rei_min)/81.4
164  IF (abs(d_rei_dt-0.71)<1.E-4) d_rei_dt = 0.71
165
166  ! Calculer l'epaisseur optique et l'emmissivite des nuages
167  ! IM inversion des DO
168
169  xflwp = 0.D0
170  xfiwp = 0.D0
171  xflwc = 0.D0
172  xfiwc = 0.D0
173
174  reliq = 0.
175  reice = 0.
176  reliq_pi = 0.
177  reice_pi = 0.
178
179  IF (iflag_t_glace.EQ.0) THEN
180    DO k = 1, klev
181      DO i = 1, klon
182        ! -layer calculation
183        rhodz(i, k) = (paprs(i,k)-paprs(i,k+1))/rg ! kg/m2
184        zrho(i, k) = pplay(i, k)/temp(i, k)/rd ! kg/m3
185        dh(i, k) = rhodz(i, k)/zrho(i, k) ! m
186        ! -Fraction of ice in cloud using a linear transition
187        icefrac_optics(i, k) = 1.0 - (temp(i,k)-t_glace_min_old)/(t_glace_max_old-t_glace_min_old)
188        icefrac_optics(i, k) = min(max(icefrac_optics(i,k),0.0), 1.0)
189        ! -IM Total Liquid/Ice water content
190        xflwc(i, k) = (1.-icefrac_optics(i,k))*radocond(i, k)
191        xfiwc(i, k) = icefrac_optics(i, k)*radocond(i, k)
192      ENDDO
193    ENDDO
194  ELSE ! of IF (iflag_t_glace.EQ.0)
195    DO k = 1, klev
196
197
198      IF (ok_new_lscp) THEN
199          CALL icefrac_lscp(klon,temp(:,k),iflag_ice_thermo,distcltop(:,k),temp_cltop(:,k),icefrac_optics(:,k),dzfice(:,k))
200      ELSE
201          CALL icefrac_lsc(klon,temp(:,k),pplay(:,k)/paprs(:,1),icefrac_optics(:,k))
202      ENDIF
203
204      DO i = 1, klon
205       
206        IF ((.NOT. ptconv(i,k)) .AND. ok_new_lscp .AND. ok_icefra_lscp) THEN
207        ! EV: take the ice fraction directly from the lscp code
208        ! consistent only for non convective grid points
209        ! critical for mixed phase clouds
210            icefrac_optics(i,k)=picefra(i,k)
211        ENDIF
212
213        ! -layer calculation
214        rhodz(i, k) = (paprs(i,k)-paprs(i,k+1))/rg ! kg/m2
215        zrho(i, k) = pplay(i, k)/temp(i, k)/rd ! kg/m3
216        dh(i, k) = rhodz(i, k)/zrho(i, k) ! m
217        ! -IM Total Liquid/Ice water content
218        xflwc(i, k) = (1.-icefrac_optics(i,k))*radocond(i, k)
219        xfiwc(i, k) = icefrac_optics(i, k)*radocond(i, k)
220      ENDDO
221    ENDDO
222  ENDIF
223
224
225
226
227
228
229  IF (ok_cdnc) THEN
230
231    ! --we compute cloud properties as a function of the aerosol load
232
233    DO k = 1, klev
234      DO i = 1, klon
235        ! Formula "D" of Boucher and Lohmann, Tellus, 1995
236        ! Cloud droplet number concentration (CDNC) is restricted
237        ! to be within [20, 1000 cm^3]
238
239        ! --pre-industrial case
240        cdnc_pi(i, k) = 10.**(bl95_b0+bl95_b1*log(max(mass_solu_aero_pi(i,k), &
241          1.E-4))/log(10.))*1.E6 !-m-3
242        cdnc_pi(i, k) = min(cdnc_max_m3, max(cdnc_min_m3,cdnc_pi(i,k)))
243
244      ENDDO
245    ENDDO
246
247    !--flag_aerosol=7 => MACv2SP climatology
248    !--in this case there is an enhancement factor
249    IF (flag_aerosol .EQ. 7) THEN
250
251      !--present-day
252      DO k = 1, klev
253        DO i = 1, klon
254          cdnc(i, k) = cdnc_pi(i,k)*dNovrN(i)
255        ENDDO
256      ENDDO
257
258    !--standard case
259    ELSE
260
261      DO k = 1, klev
262        DO i = 1, klon
263
264          ! Formula "D" of Boucher and Lohmann, Tellus, 1995
265          ! Cloud droplet number concentration (CDNC) is restricted
266          ! to be within [20, 1000 cm^3]
267
268          ! --present-day case
269          cdnc(i, k) = 10.**(bl95_b0+bl95_b1*log(max(mass_solu_aero(i,k), &
270            1.E-4))/log(10.))*1.E6 !-m-3
271          cdnc(i, k) = min(cdnc_max_m3, max(cdnc_min_m3,cdnc(i,k)))
272
273        ENDDO
274      ENDDO
275
276    ENDIF !--flag_aerosol
277
278    !--computing cloud droplet size
279    DO k = 1, klev
280      DO i = 1, klon
281
282        ! --present-day case
283        rad_chaud(i, k) = 1.1*((radocond(i,k)*pplay(i, &
284          k)/(rd*temp(i,k)))/(4./3*rpi*1000.*cdnc(i,k)))**(1./3.)
285        rad_chaud(i, k) = max(rad_chaud(i,k)*1.E6, 5.)
286
287        ! --pre-industrial case
288        rad_chaud_pi(i, k) = 1.1*((radocond(i,k)*pplay(i, &
289          k)/(rd*temp(i,k)))/(4./3.*rpi*1000.*cdnc_pi(i,k)))**(1./3.)
290        rad_chaud_pi(i, k) = max(rad_chaud_pi(i,k)*1.E6, 5.)
291
292        ! --pre-industrial case
293        ! --liquid/ice cloud water paths:
294        IF (pclc(i,k)<=seuil_neb) THEN
295
296          pcldtaupi(i, k) = 0.0
297
298        ELSE
299
300          zflwp_var = 1000.*(1.-icefrac_optics(i,k))*radocond(i, k)/pclc(i, k)* &
301            rhodz(i, k)
302          zfiwp_var = 1000.*icefrac_optics(i, k)*radocond(i, k)/pclc(i, k)*rhodz(i, k)
303          ! Calculation of ice cloud effective radius in micron
304          IF (iflag_rei .EQ. 1) THEN
305            ! when we account for precipitation in the radiation scheme,
306            ! It is recommended to use the rei formula from Sun and Rikkus 1999 with a revision
307            ! from Sun 2001 (as in the IFS model)
308            iwc=icefrac_optics(i, k)*radocond(i, k)/pclc(i,k)*zrho(i,k)*1000. !in cloud ice water content in g/m3
309            dei=(1.2351+0.0105*(temp(i,k)-273.15))*(45.8966*(iwc**0.2214) + &
310               & 0.7957*(iwc**0.2535)*(temp(i,k)-83.15))
311            !deimax=155.0
312            !deimin=20.+40*cos(abs(latitude_deg(i))/180.*RPI)
313            !Etienne: deimax and deimin controled by rei_max and rei_min in physiq.def
314            deimax=rei_max*2.0
315            deimin=2.0*rei_min+40*cos(abs(latitude_deg(i))/180.*RPI)
316            dei=min(dei,deimax)
317            dei=max(dei,deimin)
318            rei=3.*sqrt(3.)/8.*dei
319           ELSE
320            ! Default
321            ! for ice clouds: as a function of the ambiant temperature
322            ! [formula used by Iacobellis and Somerville (2000), with an
323            ! asymptotical value of 3.5 microns at T<-81.4 C added to be
324            ! consistent with observations of Heymsfield et al. 1986]:
325            ! 2011/05/24 : rei_min = 3.5 becomes a free PARAMETER as well as
326            ! rei_max=61.29
327            tc = temp(i, k) - 273.15
328            rei = d_rei_dt*tc + rei_max
329            IF (tc<=-81.4) rei = rei_min
330           ENDIF
331
332          ! -- cloud optical thickness :
333          ! [for liquid clouds, traditional formula,
334          ! for ice clouds, Ebert & Curry (1992)]
335
336          IF (zfiwp_var==0. .OR. rei<=0.) rei = 1.
337          pcldtaupi(i, k) = 3.0/2.0*zflwp_var/rad_chaud_pi(i, k) + &
338            zfiwp_var*(3.448E-03+2.431/rei)
339
340        ENDIF
341
342      ENDDO
343    ENDDO
344
345  ELSE !--not ok_cdnc
346
347    ! -prescribed cloud droplet radius
348
349    DO k = 1, min(3, klev)
350      DO i = 1, klon
351        rad_chaud(i, k) = rad_chau2
352        rad_chaud_pi(i, k) = rad_chau2
353      ENDDO
354    ENDDO
355    DO k = min(3, klev) + 1, klev
356      DO i = 1, klon
357        rad_chaud(i, k) = rad_chau1
358        rad_chaud_pi(i, k) = rad_chau1
359      ENDDO
360    ENDDO
361
362  ENDIF !--ok_cdnc
363
364  ! --computation of cloud optical depth and emissivity
365  ! --in the general case
366
367  DO k = 1, klev
368    DO i = 1, klon
369
370      IF (pclc(i,k)<=seuil_neb) THEN
371
372        ! effective cloud droplet radius (microns) for liquid water clouds:
373        ! For output diagnostics cloud droplet effective radius [um]
374        ! we multiply here with f * xl (fraction of liquid water
375        ! clouds in the grid cell) to avoid problems in the averaging of the
376        ! output.
377        ! In the output of IOIPSL, derive the REAL cloud droplet
378        ! effective radius as re/fl
379
380        fl(i, k) = seuil_neb*(1.-icefrac_optics(i,k))
381        re(i, k) = rad_chaud(i, k)*fl(i, k)
382        rel = 0.
383        rei = 0.
384        pclc(i, k) = 0.0
385        pcltau(i, k) = 0.0
386        pclemi(i, k) = 0.0
387
388      ELSE
389
390        ! -- liquid/ice cloud water paths:
391
392        zflwp_var = 1000.*(1.-icefrac_optics(i,k))*radocond(i, k)/pclc(i, k)*rhodz(i, k)
393        zfiwp_var = 1000.*icefrac_optics(i, k)*radocond(i, k)/pclc(i, k)*rhodz(i, k)
394
395        ! effective cloud droplet radius (microns) for liquid water clouds:
396        ! For output diagnostics cloud droplet effective radius [um]
397        ! we multiply here with f Effective radius of cloud droplet at top of cloud (m)* xl (fraction of liquid water
398        ! clouds in the grid cell) to avoid problems in the averaging of the
399        ! output.
400        ! In the output of IOIPSL, derive the REAL cloud droplet
401        ! effective radius as re/fl
402
403        fl(i, k) = pclc(i, k)*(1.-icefrac_optics(i,k))
404        re(i, k) = rad_chaud(i, k)*fl(i, k)
405
406        rel = rad_chaud(i, k)
407
408        ! Calculation of ice cloud effective radius in micron
409
410
411        IF (iflag_rei .GT. 0) THEN
412
413            ! when we account for precipitation in the radiation scheme,
414            ! we use the rei formula from Sun and Rikkus 1999 with a revision
415            ! from Sun 2001 (as in the IFS model)
416            iwc=icefrac_optics(i, k)*radocond(i, k)/pclc(i,k)*zrho(i,k)*1000. !in cloud ice water content in g/m3
417            dei=(1.2351+0.0105*(temp(i,k)-273.15))*(45.8966*(iwc**0.2214) + &
418               &0.7957*(iwc**0.2535)*(temp(i,k)-83.15))
419            !deimax=155.0
420            !deimin=20.+40*cos(abs(latitude_deg(i))/180.*RPI)
421            !Etienne: deimax and deimin controled by rei_max and rei_min in physiq.def
422            deimax=rei_max*2.0
423            deimin=2.0*rei_min+40*cos(abs(latitude_deg(i))/180.*RPI)
424            dei=min(dei,deimax)
425            dei=max(dei,deimin)
426            rei=3.*sqrt(3.)/8.*dei
427       
428        ELSE
429            ! Default
430            ! for ice clouds: as a function of the ambiant temperature
431            ! [formula used by Iacobellis and Somerville (2000), with an
432            ! asymptotical value of 3.5 microns at T<-81.4 C added to be
433            ! consistent with observations of Heymsfield et al. 1986]:
434            ! 2011/05/24 : rei_min = 3.5 becomes a free PARAMETER as well as
435            ! rei_max=61.29
436            tc = temp(i, k) - 273.15
437            rei = d_rei_dt*tc + rei_max
438            IF (tc<=-81.4) rei = rei_min
439        ENDIF
440        ! -- cloud optical thickness :
441        ! [for liquid clouds, traditional formula,
442        ! for ice clouds, Ebert & Curry (1992)]
443
444        IF (zflwp_var==0.) rel = 1.
445        IF (zfiwp_var==0. .OR. rei<=0.) rei = 1.
446        pcltau(i, k) = 3.0/2.0*(zflwp_var/rel) + zfiwp_var*(3.448E-03+2.431/ &
447          rei)
448
449        ! -- cloud infrared emissivity:
450        ! [the broadband infrared absorption coefficient is PARAMETERized
451        ! as a function of the effective cld droplet radius]
452        ! Ebert and Curry (1992) formula as used by Kiehl & Zender (1995):
453
454        k_ice = k_ice0 + 1.0/rei
455
456        pclemi(i, k) = 1.0 - exp(-coef_chau*zflwp_var-df*k_ice*zfiwp_var)
457
458      ENDIF
459
460      reice(i, k) = rei
461
462      xflwp(i) = xflwp(i) + xflwc(i, k)*rhodz(i, k)
463      xfiwp(i) = xfiwp(i) + xfiwc(i, k)*rhodz(i, k)
464
465    ENDDO
466  ENDDO
467
468  ! --if cloud droplet radius is fixed, then pcldtaupi=pcltau
469
470  IF (.NOT. ok_cdnc) THEN
471    DO k = 1, klev
472      DO i = 1, klon
473        pcldtaupi(i, k) = pcltau(i, k)
474        reice_pi(i, k) = reice(i, k)
475      ENDDO
476    ENDDO
477  ENDIF
478
479  DO k = 1, klev
480    DO i = 1, klon
481      reliq(i, k) = rad_chaud(i, k)
482      reliq_pi(i, k) = rad_chaud_pi(i, k)
483      reice_pi(i, k) = reice(i, k)
484    ENDDO
485  ENDDO
486
487  ! COMPUTE CLOUD LIQUID PATH AND TOTAL CLOUDINESS
488  ! IM cf. CR:test: calcul prenant ou non en compte le recouvrement
489  ! initialisations
490
491  DO i = 1, klon
492    zclear(i) = 1.
493    zcloud(i) = 0.
494    zcloudh(i) = 0.
495    zcloudm(i) = 0.
496    zcloudl(i) = 0.
497    pch(i) = 1.0
498    pcm(i) = 1.0
499    pcl(i) = 1.0
500    radocondwp(i) = 0.0
501  ENDDO
502
503  ! --calculation of liquid water path
504
505  DO k = klev, 1, -1
506    DO i = 1, klon
507      radocondwp(i) = radocondwp(i) + radocond(i, k)*rhodz(i, k)
508    ENDDO
509  ENDDO
510
511  ! --calculation of cloud properties with cloud overlap
512  ! choix de l'hypothese de recouvrement nuageuse via radopt.h (IM, 19.07.2016)
513  ! !novlp=1: max-random
514  ! !novlp=2: maximum
515  ! !novlp=3: random
516
517
518  IF (novlp==1) THEN
519    DO k = klev, 1, -1
520      DO i = 1, klon
521        zclear(i) = zclear(i)*(1.-max(pclc(i,k),zcloud(i)))/(1.-min(real( &
522          zcloud(i),kind=8),1.-zepsec))
523        pct(i) = 1. - zclear(i)
524        IF (paprs(i,k)<prmhc) THEN
525          pch(i) = pch(i)*(1.-max(pclc(i,k),zcloudh(i)))/(1.-min(real(zcloudh &
526            (i),kind=8),1.-zepsec))
527          zcloudh(i) = pclc(i, k)
528        ELSE IF (paprs(i,k)>=prmhc .AND. paprs(i,k)<prlmc) THEN
529          pcm(i) = pcm(i)*(1.-max(pclc(i,k),zcloudm(i)))/(1.-min(real(zcloudm &
530            (i),kind=8),1.-zepsec))
531          zcloudm(i) = pclc(i, k)
532        ELSE IF (paprs(i,k)>=prlmc) THEN
533          pcl(i) = pcl(i)*(1.-max(pclc(i,k),zcloudl(i)))/(1.-min(real(zcloudl &
534            (i),kind=8),1.-zepsec))
535          zcloudl(i) = pclc(i, k)
536        ENDIF
537        zcloud(i) = pclc(i, k)
538      ENDDO
539    ENDDO
540  ELSE IF (novlp==2) THEN
541    DO k = klev, 1, -1
542      DO i = 1, klon
543        zcloud(i) = max(pclc(i,k), zcloud(i))
544        pct(i) = zcloud(i)
545        IF (paprs(i,k)<prmhc) THEN
546          pch(i) = min(pclc(i,k), pch(i))
547        ELSE IF (paprs(i,k)>=prmhc .AND. paprs(i,k)<prlmc) THEN
548          pcm(i) = min(pclc(i,k), pcm(i))
549        ELSE IF (paprs(i,k)>=prlmc) THEN
550          pcl(i) = min(pclc(i,k), pcl(i))
551        ENDIF
552      ENDDO
553    ENDDO
554  ELSE IF (novlp==3) THEN
555    DO k = klev, 1, -1
556      DO i = 1, klon
557        zclear(i) = zclear(i)*(1.-pclc(i,k))
558        pct(i) = 1 - zclear(i)
559        IF (paprs(i,k)<prmhc) THEN
560          pch(i) = pch(i)*(1.0-pclc(i,k))
561        ELSE IF (paprs(i,k)>=prmhc .AND. paprs(i,k)<prlmc) THEN
562          pcm(i) = pcm(i)*(1.0-pclc(i,k))
563        ELSE IF (paprs(i,k)>=prlmc) THEN
564          pcl(i) = pcl(i)*(1.0-pclc(i,k))
565        ENDIF
566      ENDDO
567    ENDDO
568  ENDIF
569
570  DO i = 1, klon
571    pch(i) = 1. - pch(i)
572    pcm(i) = 1. - pcm(i)
573    pcl(i) = 1. - pcl(i)
574  ENDDO
575
576  ! ========================================================
577  ! DIAGNOSTICS CALCULATION FOR CMIP5 PROTOCOL
578  ! ========================================================
579  ! change by Nicolas Yan (LSCE)
580  ! Cloud Droplet Number Concentration (CDNC) : 3D variable
581  ! Fractionnal cover by liquid water cloud (LCC3D) : 3D variable
582  ! Cloud Droplet Number Concentration at top of cloud (CLDNCL) : 2D variable
583  ! Droplet effective radius at top of cloud (REFFCLWTOP) : 2D variable
584  ! Fractionnal cover by liquid water at top of clouds (LCC) : 2D variable
585
586  IF (ok_cdnc) THEN
587
588    DO k = 1, klev
589      DO i = 1, klon
590        phase3d(i, k) = 1 - icefrac_optics(i, k)
591        IF (pclc(i,k)<=seuil_neb) THEN
592          lcc3d(i, k) = seuil_neb*phase3d(i, k)
593        ELSE
594          lcc3d(i, k) = pclc(i, k)*phase3d(i, k)
595        ENDIF
596        scdnc(i, k) = lcc3d(i, k)*cdnc(i, k) ! m-3
597      ENDDO
598    ENDDO
599
600    DO i = 1, klon
601      lcc(i) = 0.
602      reffclwtop(i) = 0.
603      cldncl(i) = 0.
604      IF (novlp.EQ.3 .OR. novlp.EQ.1) tcc(i) = 1.
605      IF (novlp.EQ.2) tcc(i) = 0.
606    ENDDO
607
608    DO i = 1, klon
609      DO k = klev - 1, 1, -1 !From TOA down
610
611          ! Test, if the cloud optical depth exceeds the necessary
612          ! threshold:
613
614        IF (pcltau(i,k)>thres_tau .AND. pclc(i,k)>thres_neb) THEN
615
616          IF (novlp.EQ.2) THEN
617            IF (first) THEN
618              WRITE (*, *) 'Hypothese de recouvrement: MAXIMUM'
619              first = .FALSE.
620            ENDIF
621            flag_max = -1.
622            ftmp(i) = max(tcc(i), pclc(i,k))
623          ENDIF
624
625          IF (novlp.EQ.3) THEN
626            IF (first) THEN
627              WRITE (*, *) 'Hypothese de recouvrement: RANDOM'
628              first = .FALSE.
629            ENDIF
630            flag_max = 1.
631            ftmp(i) = tcc(i)*(1-pclc(i,k))
632          ENDIF
633
634          IF (novlp.EQ.1) THEN
635            IF (first) THEN
636              WRITE (*, *) 'Hypothese de recouvrement: MAXIMUM_ &
637                &                                             &
638                &                                          RANDOM'
639              first = .FALSE.
640            ENDIF
641            flag_max = 1.
642            ftmp(i) = tcc(i)*(1.-max(pclc(i,k),pclc(i,k+1)))/(1.-min(pclc(i, &
643              k+1),1.-thres_neb))
644          ENDIF
645          ! Effective radius of cloud droplet at top of cloud (m)
646          reffclwtop(i) = reffclwtop(i) + rad_chaud(i, k)*1.0E-06*phase3d(i, &
647            k)*(tcc(i)-ftmp(i))*flag_max
648          ! CDNC at top of cloud (m-3)
649          cldncl(i) = cldncl(i) + cdnc(i, k)*phase3d(i, k)*(tcc(i)-ftmp(i))* &
650            flag_max
651          ! Liquid Cloud Content at top of cloud
652          lcc(i) = lcc(i) + phase3d(i, k)*(tcc(i)-ftmp(i))*flag_max
653          ! Total Cloud Content at top of cloud
654          tcc(i) = ftmp(i)
655
656        ENDIF ! is there a visible, not-too-small cloud?
657      ENDDO ! loop over k
658
659      IF (novlp.EQ.3 .OR. novlp.EQ.1) tcc(i) = 1. - tcc(i)
660
661    ENDDO ! loop over i
662
663    ! ! Convective and Stratiform Cloud Droplet Effective Radius (REFFCLWC
664    ! REFFCLWS)
665    DO i = 1, klon
666      DO k = 1, klev
667        ! Weight to be used for outputs: eau_liquide*couverture nuageuse
668        lcc3dcon(i, k) = rnebcon(i, k)*phase3d(i, k)*ccwcon(i, k) ! eau liquide convective
669        lcc3dstra(i, k) = pclc(i, k)*radocond(i, k)*phase3d(i, k)
670        lcc3dstra(i, k) = lcc3dstra(i, k) - lcc3dcon(i, k) ! eau liquide stratiforme
671        lcc3dstra(i, k) = max(lcc3dstra(i,k), 0.0)
672        !FC pour la glace (CAUSES)
673        icc3dcon(i, k) = rnebcon(i, k)*(1-phase3d(i, k))*ccwcon(i, k) !  glace convective
674        icc3dstra(i, k)= pclc(i, k)*radocond(i, k)*(1-phase3d(i, k))
675        icc3dstra(i, k) = icc3dstra(i, k) - icc3dcon(i, k) ! glace stratiforme
676        icc3dstra(i, k) = max( icc3dstra(i, k), 0.0)
677        !FC (CAUSES)
678
679        ! Compute cloud droplet radius as above in meter
680        radius = 1.1*((radocond(i,k)*pplay(i,k)/(rd*temp(i,k)))/(4./3*rpi*1000.* &
681          cdnc(i,k)))**(1./3.)
682        radius = max(radius, 5.E-6)
683        ! Convective Cloud Droplet Effective Radius (REFFCLWC) : variable 3D
684        reffclwc(i, k) = radius
685        reffclwc(i, k) = reffclwc(i, k)*lcc3dcon(i, k)
686        ! Stratiform Cloud Droplet Effective Radius (REFFCLWS) : variable 3D
687        reffclws(i, k) = radius
688        reffclws(i, k) = reffclws(i, k)*lcc3dstra(i, k)
689      ENDDO !klev
690    ENDDO !klon
691
692    ! Column Integrated Cloud Droplet Number (CLDNVI) : variable 2D
693
694    DO i = 1, klon
695      cldnvi(i) = 0.
696      lcc_integrat(i) = 0.
697      height(i) = 0.
698      DO k = 1, klev
699        cldnvi(i) = cldnvi(i) + cdnc(i, k)*lcc3d(i, k)*dh(i, k)
700        lcc_integrat(i) = lcc_integrat(i) + lcc3d(i, k)*dh(i, k)
701        height(i) = height(i) + dh(i, k)
702      ENDDO ! klev
703      lcc_integrat(i) = lcc_integrat(i)/height(i)
704      IF (lcc_integrat(i)<=1.0E-03) THEN
705        cldnvi(i) = cldnvi(i)*lcc(i)/seuil_neb
706      ELSE
707        cldnvi(i) = cldnvi(i)*lcc(i)/lcc_integrat(i)
708      ENDIF
709    ENDDO ! klon
710
711    DO i = 1, klon
712      DO k = 1, klev
713        IF (scdnc(i,k)<=0.0) scdnc(i, k) = 0.0
714        IF (reffclws(i,k)<=0.0) reffclws(i, k) = 0.0
715        IF (reffclwc(i,k)<=0.0) reffclwc(i, k) = 0.0
716        IF (lcc3d(i,k)<=0.0) lcc3d(i, k) = 0.0
717        IF (lcc3dcon(i,k)<=0.0) lcc3dcon(i, k) = 0.0
718        IF (lcc3dstra(i,k)<=0.0) lcc3dstra(i, k) = 0.0
719!FC (CAUSES)
720        IF (icc3dcon(i,k)<=0.0) icc3dcon(i, k) = 0.0
721        IF (icc3dstra(i,k)<=0.0) icc3dstra(i, k) = 0.0
722!FC (CAUSES)
723      ENDDO
724      IF (reffclwtop(i)<=0.0) reffclwtop(i) = 0.0
725      IF (cldncl(i)<=0.0) cldncl(i) = 0.0
726      IF (cldnvi(i)<=0.0) cldnvi(i) = 0.0
727      IF (lcc(i)<=0.0) lcc(i) = 0.0
728    ENDDO
729
730  ENDIF !ok_cdnc
731
732  first=.false. !to be sure
733
734  RETURN
735
736END SUBROUTINE cloud_optics_prop
737
738END MODULE lmdz_cloud_optics_prop
Note: See TracBrowser for help on using the repository browser.