source: LMDZ6/branches/LMDZ-COSP/libf/phylmd/lmdz_cloud_optics_prop.f90 @ 5926

Last change on this file since 5926 was 5828, checked in by rkazeroni, 3 months ago

For GPU porting of call_cloud_optics_prop routine:

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