source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_cloud.F90 @ 5408

Last change on this file since 5408 was 4911, checked in by idelkadi, 8 months ago

Introduction of 2 routines for LMDZ-ECRAD :

  • One for LMDZ (calcul_cloud_overlap_decorr_len.F90) to calculate the decorrelation length (Ld) in the case of "Exp-Ran" cloud overlap.
    1. Ld=constant
    2. Ld=linear function of the absolute value of latitude (Shonk and Hogan 2010)
    3. Ld=function of vertical level
  • One for Ecrad (set_overlap_param_var2D in radiation_cloud.F90) : Compute and store the overlap parameter from the provided overlap decorrelation length, which depends on vertical level.

Translated with DeepL.com (free version)

File size: 34.7 KB
Line 
1! radiation_cloud.F90 - Derived type to store cloud/precip properties
2!
3! (C) Copyright 2014- ECMWF.
4!
5! This software is licensed under the terms of the Apache Licence Version 2.0
6! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
7!
8! In applying this licence, ECMWF does not waive the privileges and immunities
9! granted to it by virtue of its status as an intergovernmental organisation
10! nor does it submit to any jurisdiction.
11!
12! Author:  Robin Hogan
13! Email:   r.j.hogan@ecmwf.int
14!
15! Modifications
16!   2019-01-14  R. Hogan  Added inv_inhom_effective_size variable
17!   2019-01-14  R. Hogan  Added out_of_physical_bounds routine
18!   2019-06-14  R. Hogan  Added capability to store any number of cloud/precip types
19
20module radiation_cloud
21
22  use parkind1, only : jprb
23
24  implicit none
25  public
26
27  !---------------------------------------------------------------------
28  ! The intention is that all variables describing clouds and
29  ! radiatively-active precipitation are contained in this derived
30  ! type, and if cloud variables are to be added in future, they can
31  ! be added to this type without requiring extra variables to be
32  ! passed between subroutines elsewhere in the program.
33  type cloud_type
34    ! For maximum flexibility, an arbitrary number "ntype" of
35    ! hydrometeor types can be stored, dimensioned (ncol,nlev,ntype)
36    integer                                   :: ntype = 0
37    real(jprb), allocatable, dimension(:,:,:) :: &
38         &  mixing_ratio, &  ! mass mixing ratio (kg/kg)
39         &  effective_radius ! (m)
40
41    ! For backwards compatibility, we also allow for the two
42    ! traditional cloud types, liquid cloud droplets and ice cloud
43    ! particles, dimensioned (ncol,nlev)
44    real(jprb), pointer, dimension(:,:) :: &
45         &  q_liq,  q_ice,  & ! mass mixing ratio (kg/kg)
46         &  re_liq, re_ice    ! effective radius (m)
47
48    ! For the moment, the different types of hydrometeor are assumed
49    ! to be mixed with each other, so there is just one cloud fraction
50    ! variable varying from 0 to 1
51    real(jprb), allocatable, dimension(:,:) :: fraction
52
53    ! The fractional standard deviation of cloud optical depth in the
54    ! cloudy part of the gridbox.  In the Tripleclouds representation
55    ! of cloud inhomogeneity, this is implemented by splitting the
56    ! cloudy part of the gridbox into two equal-area regions, one
57    ! with the cloud optical depth scaled by 1+fractional_std and the
58    ! other scaled by 1-fractional_std. This variable is dimensioned
59    ! (ncol,nlev)
60    real(jprb), allocatable, dimension(:,:) :: fractional_std
61
62    ! The inverse of the effective horizontal size of the clouds in
63    ! the gridbox, used to compute the cloud edge length per unit
64    ! gridbox area for use in representing 3D effects. This variable
65    ! is dimensioned (ncol,nlev).
66    real(jprb), allocatable, dimension(:,:) :: inv_cloud_effective_size ! m-1
67
68    ! Similarly for the in-cloud heterogeneities, used to compute the
69    ! edge length between the optically thin and thick cloudy regions
70    ! of the gridbox.
71    real(jprb), allocatable, dimension(:,:) :: inv_inhom_effective_size ! m-1
72
73    ! The following variable describes the overlap of cloud boundaries
74    ! in adjacent layers, with dimensions (ncol,nlev-1): 1 corresponds
75    ! to maximum overlap and 0 to random overlap. Depending on the
76    ! ecRad configuration, it may be the "alpha" overlap parameter of
77    ! Hogan and Illingworth (2000) or the "beta" overlap parameter of
78    ! Shonk et al. (2010).
79    real(jprb), allocatable, dimension(:,:) :: overlap_param
80
81  contains
82    procedure :: allocate   => allocate_cloud_arrays
83    procedure :: deallocate => deallocate_cloud_arrays
84    procedure :: set_overlap_param_fix
85    procedure :: set_overlap_param_var
86    generic   :: set_overlap_param => set_overlap_param_fix, set_overlap_param_var
87    procedure :: set_overlap_param_var2D
88    procedure :: set_overlap_param_approx
89    procedure :: create_fractional_std
90    procedure :: create_inv_cloud_effective_size
91    procedure :: create_inv_cloud_effective_size_eta
92    procedure :: param_cloud_effective_separation_eta
93    procedure :: crop_cloud_fraction
94    procedure :: out_of_physical_bounds
95
96  end type cloud_type
97
98contains
99
100  !---------------------------------------------------------------------
101  ! Allocate arrays for describing clouds and precipitation, although
102  ! in the offline code these are allocated when they are read from
103  ! the NetCDF file
104  subroutine allocate_cloud_arrays(this, ncol, nlev, ntype, use_inhom_effective_size)
105
106    use yomhook,     only : lhook, dr_hook, jphook
107
108    class(cloud_type), intent(inout), target :: this
109    integer, intent(in)              :: ncol   ! Number of columns
110    integer, intent(in)              :: nlev   ! Number of levels
111    ! Number of cloud/precip particle types.  If not present then the
112    ! older cloud behaviour is assumed: two types are present, (1)
113    ! liquid and (2) ice, and they can be accessed via q_liq, q_ice,
114    ! re_liq and re_ice.
115    integer, intent(in), optional    :: ntype
116    logical, intent(in), optional    :: use_inhom_effective_size
117
118    real(jphook) :: hook_handle
119
120    if (lhook) call dr_hook('radiation_cloud:allocate',0,hook_handle)
121
122    if (present(ntype)) then
123      this%ntype = ntype
124    else
125      this%ntype = 2
126    end if
127    allocate(this%mixing_ratio(ncol,nlev,this%ntype))
128    allocate(this%effective_radius(ncol,nlev,this%ntype))
129    nullify(this%q_liq)
130    nullify(this%q_ice)
131    nullify(this%re_liq)
132    nullify(this%re_ice)
133    if (.not. present(ntype)) then
134      ! Older interface in which only liquid and ice are supported
135      this%q_liq  => this%mixing_ratio(:,:,1)
136      this%q_ice  => this%mixing_ratio(:,:,2)
137      this%re_liq => this%effective_radius(:,:,1)
138      this%re_ice => this%effective_radius(:,:,2)
139    end if
140
141    allocate(this%fraction(ncol,nlev))
142    allocate(this%overlap_param(ncol,nlev-1))
143    allocate(this%fractional_std(ncol,nlev))
144    allocate(this%inv_cloud_effective_size(ncol,nlev))
145
146    if (present(use_inhom_effective_size)) then
147      if (use_inhom_effective_size) then
148        allocate(this%inv_inhom_effective_size(ncol,nlev))
149      end if
150    end if
151
152    if (lhook) call dr_hook('radiation_cloud:allocate',1,hook_handle)
153
154  end subroutine allocate_cloud_arrays
155
156
157  !---------------------------------------------------------------------
158  ! Deallocate arrays
159  subroutine deallocate_cloud_arrays(this)
160
161    use yomhook,     only : lhook, dr_hook, jphook
162
163    class(cloud_type), intent(inout) :: this
164
165    real(jphook) :: hook_handle
166
167    if (lhook) call dr_hook('radiation_cloud:deallocate',0,hook_handle)
168
169    nullify(this%q_liq)
170    nullify(this%q_ice)
171    nullify(this%re_liq)
172    nullify(this%re_ice)
173
174    if (allocated(this%mixing_ratio))     deallocate(this%mixing_ratio)
175    if (allocated(this%effective_radius)) deallocate(this%effective_radius)
176    if (allocated(this%fraction))         deallocate(this%fraction)
177    if (allocated(this%overlap_param))    deallocate(this%overlap_param)
178    if (allocated(this%fractional_std))   deallocate(this%fractional_std)
179    if (allocated(this%inv_cloud_effective_size)) &
180         &  deallocate(this%inv_cloud_effective_size)
181    if (allocated(this%inv_inhom_effective_size)) &
182         &  deallocate(this%inv_inhom_effective_size)
183
184    if (lhook) call dr_hook('radiation_cloud:deallocate',1,hook_handle)
185
186  end subroutine deallocate_cloud_arrays
187
188
189  !---------------------------------------------------------------------
190  ! Compute and store the overlap parameter from the provided overlap
191  ! decorrelation length (in metres).  If istartcol and/or iendcol are
192  ! provided then only columns in this range are computed.  If the
193  ! overlap_param array has not been allocated then it will be
194  ! allocated to be of the correct size relative to the pressure
195  ! field. This version assumes a fixed decorrelation_length for all
196  ! columns.
197  subroutine set_overlap_param_fix(this, thermodynamics, decorrelation_length, &
198       &  istartcol, iendcol)
199
200    use yomhook,                  only : lhook, dr_hook, jphook
201    use radiation_thermodynamics, only : thermodynamics_type
202    use radiation_constants,      only : GasConstantDryAir, AccelDueToGravity
203
204    class(cloud_type),         intent(inout) :: this
205    type(thermodynamics_type), intent(in)    :: thermodynamics
206    real(jprb),                intent(in)    :: decorrelation_length ! m
207    integer,         optional, intent(in)    :: istartcol, iendcol
208
209    ! Ratio of gas constant for dry air to acceleration due to gravity
210    real(jprb), parameter :: R_over_g = GasConstantDryAir / AccelDueToGravity
211
212    ! Process only columns i1 to i2, which will be istartcol to
213    ! iendcol if they were provided
214    integer :: i1, i2
215
216    integer :: ncol, nlev
217
218    integer :: jcol, jlev
219
220    real(jphook) :: hook_handle
221
222    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_fix',0,hook_handle)
223
224    ! Pressure at half-levels, pressure_hl, is defined at nlev+1
225    ! points
226    ncol = size(thermodynamics%pressure_hl,dim=1)
227    nlev = size(thermodynamics%pressure_hl,dim=2)-1
228
229    if (present(istartcol)) then
230      i1 = istartcol
231    else
232      i1 = 1
233    end if
234
235    if (present(iendcol)) then
236      i2 = iendcol
237    else
238      i2 = ncol
239    end if
240
241    if (.not. allocated(this%overlap_param)) then
242      ! If pressure is of size (ncol,nlev+1) then overlap_param is of
243      ! size (ncol,nlev-1), since overlap parameter is only defined here
244      ! for interfaces between model layers, not for the interface to
245      ! space or the surface
246      allocate(this%overlap_param(ncol, nlev-1))
247    end if
248
249    if (thermodynamics%pressure_hl(i1,2) > thermodynamics%pressure_hl(i1,1)) then
250      ! Pressure is increasing with index (order of layers is
251      ! top-of-atmosphere to surface). In case pressure_hl(:,1)=0, we
252      ! don't take the logarithm of the first pressure in each column.
253      do jcol = i1,i2
254        this%overlap_param(jcol,1) = exp(-(R_over_g/decorrelation_length) &
255             &                            * thermodynamics%temperature_hl(jcol,2) &
256             &                            *log(thermodynamics%pressure_hl(jcol,3) &
257             &                                /thermodynamics%pressure_hl(jcol,2)))
258      end do
259
260      do jlev = 2,nlev-1
261        do jcol = i1,i2
262          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length) &
263              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
264              &                            *log(thermodynamics%pressure_hl(jcol,jlev+2) &
265              &                                /thermodynamics%pressure_hl(jcol,jlev)))
266        end do
267      end do
268
269    else
270       ! Pressure is decreasing with index (order of layers is surface
271       ! to top-of-atmosphere).  In case pressure_hl(:,nlev+1)=0, we
272       ! don't take the logarithm of the last pressure in each column.
273      do jlev = 1,nlev-2
274        do jcol = i1,i2
275          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length) &
276              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
277              &                            *log(thermodynamics%pressure_hl(jcol,jlev) &
278              &                                /thermodynamics%pressure_hl(jcol,jlev+2)))
279        end do
280      end do
281
282      do jcol = i1,i2
283        this%overlap_param(jcol,nlev-1) = exp(-(R_over_g/decorrelation_length) &
284            &                            * thermodynamics%temperature_hl(jcol,nlev) &
285            &                            *log(thermodynamics%pressure_hl(jcol,nlev-1) &
286            &                                /thermodynamics%pressure_hl(jcol,nlev)))
287      end do
288    end if
289
290    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_fix',1,hook_handle)
291
292  end subroutine set_overlap_param_fix
293
294
295  !---------------------------------------------------------------------
296  ! Compute and store the overlap parameter from the provided overlap
297  ! decorrelation length (in metres), which may vary with column. Only
298  ! columns from istartcol to iendcol are computed.  If the
299  ! overlap_param array has not been allocated then it will be
300  ! allocated to be of the correct size relative to the pressure
301  ! field.
302  subroutine set_overlap_param_var(this, thermodynamics, decorrelation_length, &
303       &                           istartcol, iendcol)
304
305    use yomhook,                  only : lhook, dr_hook, jphook
306    use radiation_thermodynamics, only : thermodynamics_type
307    use radiation_constants,      only : GasConstantDryAir, AccelDueToGravity
308
309    class(cloud_type),         intent(inout) :: this
310    type(thermodynamics_type), intent(in)    :: thermodynamics
311    integer,                   intent(in)    :: istartcol, iendcol
312    real(jprb),                intent(in)    :: decorrelation_length(istartcol:iendcol) ! m
313
314    ! Ratio of gas constant for dry air to acceleration due to gravity
315    real(jprb), parameter :: R_over_g = GasConstantDryAir / AccelDueToGravity
316
317    integer :: ncol, nlev
318
319    integer :: jcol, jlev
320
321    real(jphook) :: hook_handle
322
323    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_var',0,hook_handle)
324
325    ! Pressure at half-levels, pressure_hl, is defined at nlev+1
326    ! points
327    ncol = size(thermodynamics%pressure_hl,dim=1)
328    nlev = size(thermodynamics%pressure_hl,dim=2)-1
329
330    if (.not. allocated(this%overlap_param)) then
331      ! If pressure is of size (ncol,nlev+1) then overlap_param is of
332      ! size (ncol,nlev-1), since overlap parameter is only defined here
333      ! for interfaces between model layers, not for the interface to
334      ! space or the surface
335      allocate(this%overlap_param(ncol, nlev-1))
336    end if
337
338    if (thermodynamics%pressure_hl(istartcol,2) > thermodynamics%pressure_hl(istartcol,1)) then
339      ! Pressure is increasing with index (order of layers is
340      ! top-of-atmosphere to surface). In case pressure_hl(:,1)=0, we
341      ! don't take the logarithm of the first pressure in each column.
342      do jcol = istartcol,iendcol
343        this%overlap_param(jcol,1) = exp(-(R_over_g/decorrelation_length(jcol)) &
344             &                            * thermodynamics%temperature_hl(jcol,2) &
345             &                            *log(thermodynamics%pressure_hl(jcol,3) &
346             &                                /thermodynamics%pressure_hl(jcol,2)))
347      end do
348
349      do jlev = 2,nlev-1
350        do jcol = istartcol,iendcol
351          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length(jcol)) &
352              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
353              &                            *log(thermodynamics%pressure_hl(jcol,jlev+2) &
354              &                                /thermodynamics%pressure_hl(jcol,jlev)))
355        end do
356      end do
357
358    else
359       ! Pressure is decreasing with index (order of layers is surface
360       ! to top-of-atmosphere).  In case pressure_hl(:,nlev+1)=0, we
361       ! don't take the logarithm of the last pressure in each column.
362      do jlev = 1,nlev-2
363        do jcol = istartcol,iendcol
364          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length(jcol)) &
365              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
366              &                            *log(thermodynamics%pressure_hl(jcol,jlev) &
367              &                                /thermodynamics%pressure_hl(jcol,jlev+2)))
368        end do
369      end do
370
371      do jcol = istartcol,iendcol
372        this%overlap_param(jcol,nlev-1) = exp(-(R_over_g/decorrelation_length(jcol)) &
373            &                            * thermodynamics%temperature_hl(jcol,nlev) &
374            &                            *log(thermodynamics%pressure_hl(jcol,nlev-1) &
375            &                                /thermodynamics%pressure_hl(jcol,nlev)))
376      end do
377    end if
378
379    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_var',1,hook_handle)
380
381  end subroutine set_overlap_param_var
382
383! AI 04 2024 variation de la longueur Ld en fonction de la verticale 
384  subroutine set_overlap_param_var2D(this, thermodynamics, decorrelation_length, &
385       &                           klev, istartcol, iendcol)
386
387    use yomhook,                  only : lhook, dr_hook, jphook
388    use radiation_thermodynamics, only : thermodynamics_type
389    use radiation_constants,      only : GasConstantDryAir, AccelDueToGravity
390
391    integer,                   intent(in)    :: klev
392    class(cloud_type),         intent(inout) :: this
393    type(thermodynamics_type), intent(in)    :: thermodynamics
394    integer,                   intent(in)    :: istartcol, iendcol
395    real(jprb),                intent(in)    :: decorrelation_length(istartcol:iendcol,klev) ! m
396
397    ! Ratio of gas constant for dry air to acceleration due to gravity
398    real(jprb), parameter :: R_over_g = GasConstantDryAir / AccelDueToGravity
399
400    integer :: ncol, nlev
401
402    integer :: jcol, jlev
403
404    real(jphook) :: hook_handle
405
406    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_var',0,hook_handle)
407
408    ! Pressure at half-levels, pressure_hl, is defined at nlev+1
409    ! points
410    ncol = size(thermodynamics%pressure_hl,dim=1)
411    nlev = size(thermodynamics%pressure_hl,dim=2)-1
412
413    if (.not. allocated(this%overlap_param)) then
414      ! If pressure is of size (ncol,nlev+1) then overlap_param is of
415      ! size (ncol,nlev-1), since overlap parameter is only defined here
416      ! for interfaces between model layers, not for the interface to
417      ! space or the surface
418      allocate(this%overlap_param(ncol, nlev-1))
419    end if
420
421    if (thermodynamics%pressure_hl(istartcol,2) > thermodynamics%pressure_hl(istartcol,1)) then
422      ! Pressure is increasing with index (order of layers is
423      ! top-of-atmosphere to surface). In case pressure_hl(:,1)=0, we
424      ! don't take the logarithm of the first pressure in each column.
425      do jcol = istartcol,iendcol
426        this%overlap_param(jcol,1) = exp(-(R_over_g/decorrelation_length(jcol,1)) &
427             &                            * thermodynamics%temperature_hl(jcol,2) &
428             &                            *log(thermodynamics%pressure_hl(jcol,3) &
429             &                                /thermodynamics%pressure_hl(jcol,2)))
430      end do
431
432      do jlev = 2,nlev-1
433        do jcol = istartcol,iendcol
434          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length(jcol,jlev)) &
435              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
436              &                            *log(thermodynamics%pressure_hl(jcol,jlev+2) &
437              &                                /thermodynamics%pressure_hl(jcol,jlev)))
438        end do
439      end do
440
441    else
442       ! Pressure is decreasing with index (order of layers is surface
443       ! to top-of-atmosphere).  In case pressure_hl(:,nlev+1)=0, we
444       ! don't take the logarithm of the last pressure in each column.
445      do jlev = 1,nlev-2
446        do jcol = istartcol,iendcol
447          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length(jcol,jlev)) &
448              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
449              &                            *log(thermodynamics%pressure_hl(jcol,jlev) &
450              &                                /thermodynamics%pressure_hl(jcol,jlev+2)))
451        end do
452      end do
453
454      do jcol = istartcol,iendcol
455      ! AI ATTENTION a verifier decorrelation_length(jcol,nlev-1) ou nlev
456        this%overlap_param(jcol,nlev-1) = exp(-(R_over_g/decorrelation_length(jcol,nlev-1)) &
457            &                            * thermodynamics%temperature_hl(jcol,nlev) &
458            &                            *log(thermodynamics%pressure_hl(jcol,nlev-1) &
459            &                                /thermodynamics%pressure_hl(jcol,nlev)))
460      end do
461    end if
462
463    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_var',1,hook_handle)
464
465  end subroutine set_overlap_param_var2D
466
467  !---------------------------------------------------------------------
468  ! Compute and store the overlap parameter from the provided overlap
469  ! decorrelation length (in metres).  If istartcol and/or iendcol are
470  ! provided then only columns in this range are computed.  If the
471  ! overlap_param array has not been allocated then it will be
472  ! allocated to be of the correct size relative to the pressure
473  ! field. This is the APPROXIMATE method as it assumes a fixed
474  ! atmospheric scale height, which leads to differences particularly
475  ! in low cloud.
476  subroutine set_overlap_param_approx(this, thermodynamics, decorrelation_length, &
477       &  istartcol, iendcol)
478
479    use yomhook,                  only : lhook, dr_hook, jphook
480    use radiation_thermodynamics, only : thermodynamics_type
481
482    class(cloud_type),         intent(inout) :: this
483    type(thermodynamics_type), intent(in)    :: thermodynamics
484    real(jprb),                intent(in)    :: decorrelation_length ! m
485    integer,         optional, intent(in)    :: istartcol, iendcol
486
487    ! To convert decorrelation length (m) to overlap parameter between
488    ! layers, we need an estimate for the thickness of the layer. This
489    ! is found using the pressure difference between the edges of the
490    ! layer, along with the approximate scale height of the atmosphere
491    ! (m) given here:
492    real(jprb), parameter :: scale_height = 8000.0_jprb
493
494    ! Process only columns i1 to i2, which will be istartcol to
495    ! iendcol if they were provided
496    integer :: i1, i2
497
498    integer :: ncol, nlev
499
500    real(jphook) :: hook_handle
501
502    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_approx',0,hook_handle)
503
504    ! Pressure at half-levels, pressure_hl, is defined at nlev+1
505    ! points
506    ncol = size(thermodynamics%pressure_hl,dim=1)
507    nlev = size(thermodynamics%pressure_hl,dim=2)-1
508
509    if (present(istartcol)) then
510      i1 = istartcol
511    else
512      i1 = 1
513    end if
514
515    if (present(iendcol)) then
516      i2 = iendcol
517    else
518      i2 = ncol
519    end if
520
521    if (.not. allocated(this%overlap_param)) then
522      ! If pressure is of size (ncol,nlev+1) then overlap_param is of
523      ! size (ncol,nlev-1), since overlap parameter is only defined here
524      ! for interfaces between model layers, not for the interface to
525      ! space or the surface
526      allocate(this%overlap_param(ncol, nlev-1))
527    end if
528
529    if (thermodynamics%pressure_hl(i1,2) > thermodynamics%pressure_hl(i1,1)) then
530       ! Pressure is increasing with index (order of layers is
531       ! top-of-atmosphere to surface). In case pressure_hl(:,1)=0, we
532       ! don't take the logarithm of the first pressure in each
533       ! column.
534       this%overlap_param(i1:i2,:) = exp(-(scale_height/decorrelation_length) &
535            &  * ( log(thermodynamics%pressure_hl(i1:i2,3:nlev+1) &
536            &         /thermodynamics%pressure_hl(i1:i2,2:nlev  )) ) )
537    else
538       ! Pressure is decreasing with index (order of layers is surface
539       ! to top-of-atmosphere).  In case pressure_hl(:,nlev+1)=0, we
540       ! don't take the logarithm of the last pressure in each column.
541       this%overlap_param(i1:i2,:) = exp(-(scale_height/decorrelation_length) &
542            &  * ( log(thermodynamics%pressure_hl(i1:i2,1:nlev-1) &
543            &         /thermodynamics%pressure_hl(i1:i2,2:nlev  )) ) )
544    end if
545
546    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_approx',1,hook_handle)
547
548  end subroutine set_overlap_param_approx
549
550
551  !---------------------------------------------------------------------
552  ! Create a matrix of constant fractional standard deviations
553  ! (dimensionless)
554  subroutine create_fractional_std(this, ncol, nlev, frac_std)
555
556    use yomhook,                  only : lhook, dr_hook, jphook
557
558    class(cloud_type), intent(inout) :: this
559    integer,           intent(in)    :: ncol, nlev
560    real(jprb),        intent(in)    :: frac_std
561
562    real(jphook) :: hook_handle
563
564    if (lhook) call dr_hook('radiation_cloud:create_fractional_std',0,hook_handle)
565
566    if (allocated(this%fractional_std)) then
567       deallocate(this%fractional_std)
568    end if
569   
570    allocate(this%fractional_std(ncol, nlev))
571
572    this%fractional_std = frac_std
573
574    if (lhook) call dr_hook('radiation_cloud:create_fractional_std',1,hook_handle)
575
576  end subroutine create_fractional_std
577
578
579  !---------------------------------------------------------------------
580  ! Create a matrix of constant inverse cloud effective size (m-1)
581  subroutine create_inv_cloud_effective_size(this, ncol, nlev, inv_eff_size)
582
583    use yomhook,                  only : lhook, dr_hook, jphook
584
585    class(cloud_type), intent(inout) :: this
586    integer,           intent(in)    :: ncol, nlev
587    real(jprb),        intent(in)    :: inv_eff_size
588
589    real(jphook) :: hook_handle
590
591    if (lhook) call dr_hook('radiation_cloud:create_inv_cloud_effective_size',0,hook_handle)
592
593    if (allocated(this%inv_cloud_effective_size)) then
594       deallocate(this%inv_cloud_effective_size)
595    end if
596   
597    allocate(this%inv_cloud_effective_size(ncol, nlev))
598
599    this%inv_cloud_effective_size = inv_eff_size
600
601    if (lhook) call dr_hook('radiation_cloud:create_inv_cloud_effective_size',1,hook_handle)
602
603  end subroutine create_inv_cloud_effective_size
604
605
606  !---------------------------------------------------------------------
607  ! Create a matrix of inverse cloud effective size (m-1) according to
608  ! the value of eta (=pressure divided by surface pressure)
609  subroutine create_inv_cloud_effective_size_eta(this, ncol, nlev, &
610       &  pressure_hl, inv_eff_size_low, inv_eff_size_mid, inv_eff_size_high, &
611       &  eta_low_mid, eta_mid_high, istartcol, iendcol)
612
613    use yomhook,                  only : lhook, dr_hook, jphook
614
615    class(cloud_type), intent(inout) :: this
616    integer,           intent(in)    :: ncol, nlev
617    ! Pressure on half levels (Pa)
618    real(jprb),        intent(in)    :: pressure_hl(:,:)
619    ! Inverse effective size for low, mid and high cloud (m-1)
620    real(jprb),        intent(in)    :: inv_eff_size_low
621    real(jprb),        intent(in)    :: inv_eff_size_mid
622    real(jprb),        intent(in)    :: inv_eff_size_high
623    ! Eta values at low-mid and mid-high interfaces
624    real(jprb),        intent(in)    :: eta_low_mid, eta_mid_high
625    integer, optional, intent(in)    :: istartcol, iendcol
626
627    ! Ratio of layer midpoint pressure to surface pressure
628    real(jprb) :: eta(nlev)
629
630    ! Indices of column, level and surface half-level
631    integer :: jcol, isurf
632
633    ! Local values of istartcol, iendcol
634    integer :: i1, i2
635
636    real(jphook) :: hook_handle
637
638    if (lhook) call dr_hook('radiation_cloud:create_inv_cloud_effective_size_eta',0,hook_handle)
639
640    if (allocated(this%inv_cloud_effective_size)) then
641      deallocate(this%inv_cloud_effective_size)
642    end if
643   
644    allocate(this%inv_cloud_effective_size(ncol, nlev))
645
646    if (present(istartcol)) then
647      i1 = istartcol
648    else
649      i1 = 1
650    end if
651
652    if (present(iendcol)) then
653      i2 = iendcol
654    else
655      i2 = ncol
656    end if
657
658    ! Locate the surface half-level
659    if (pressure_hl(1,1) > pressure_hl(1,2)) then
660      isurf = 1
661    else
662      isurf = nlev+1
663    end if
664
665    do jcol = i1,i2
666      eta = (pressure_hl(jcol,1:nlev)+pressure_hl(jcol,2:nlev+1)) &
667           &  * (0.5_jprb / pressure_hl(jcol,isurf))
668      where (eta > eta_low_mid)
669        this%inv_cloud_effective_size(jcol,:) = inv_eff_size_low
670      elsewhere (eta > eta_mid_high)
671        this%inv_cloud_effective_size(jcol,:) = inv_eff_size_mid
672      elsewhere
673        this%inv_cloud_effective_size(jcol,:) = inv_eff_size_high
674      end where
675    end do
676
677    if (lhook) call dr_hook('radiation_cloud:create_inv_cloud_effective_size_eta',1,hook_handle)
678
679  end subroutine create_inv_cloud_effective_size_eta
680
681
682  !---------------------------------------------------------------------
683  ! Create a matrix of inverse cloud and inhomogeneity effective size
684  ! (m-1) parameterized according to the value of eta (=pressure
685  ! divided by surface pressure): effective_separation =
686  ! coeff_a + coeff_b*exp(-(eta**power)). 
687  subroutine param_cloud_effective_separation_eta(this, ncol, nlev, &
688       &  pressure_hl, separation_surf, separation_toa, power, &
689       &  inhom_separation_factor, istartcol, iendcol)
690
691    use yomhook,                  only : lhook, dr_hook, jphook
692
693    class(cloud_type), intent(inout) :: this
694    integer,           intent(in)    :: ncol, nlev
695    ! Pressure on half levels (Pa)
696    real(jprb),        intent(in)    :: pressure_hl(:,:)
697    ! Separation distances at surface and top-of-atmosphere, and power
698    ! on eta
699    real(jprb),           intent(in) :: separation_surf ! m
700    real(jprb),           intent(in) :: separation_toa ! m
701    real(jprb),           intent(in) :: power
702    real(jprb), optional, intent(in) :: inhom_separation_factor
703    integer,    optional, intent(in) :: istartcol, iendcol
704
705    ! Ratio of layer midpoint pressure to surface pressure
706    real(jprb) :: eta(nlev)
707
708    ! Effective cloud separation (m)
709    real(jprb) :: eff_separation(nlev)
710
711    ! Coefficients used to compute effective separation distance
712    real(jprb) :: coeff_e, coeff_a, coeff_b, inhom_sep_factor
713
714    ! Indices of column, level and surface half-level
715    integer :: jcol, isurf
716
717    ! Local values of istartcol, iendcol
718    integer :: i1, i2
719
720    real(jphook) :: hook_handle
721
722    if (lhook) call dr_hook('radiation_cloud:param_cloud_effective_separation_eta',0,hook_handle)
723
724    if (present(inhom_separation_factor)) then
725      inhom_sep_factor = inhom_separation_factor
726    else
727      inhom_sep_factor = 1.0_jprb
728    end if
729
730    coeff_e = 1.0_jprb - exp(-1.0_jprb)
731    coeff_b = (separation_toa - separation_surf) / coeff_e
732    coeff_a = separation_toa - coeff_b
733
734    if (allocated(this%inv_cloud_effective_size)) then
735      deallocate(this%inv_cloud_effective_size)
736    end if
737     if (allocated(this%inv_inhom_effective_size)) then
738      deallocate(this%inv_inhom_effective_size)
739    end if
740   
741    allocate(this%inv_cloud_effective_size(ncol, nlev))
742    allocate(this%inv_inhom_effective_size(ncol, nlev))
743
744    if (present(istartcol)) then
745      i1 = istartcol
746    else
747      i1 = 1
748    end if
749
750    if (present(iendcol)) then
751      i2 = iendcol
752    else
753      i2 = ncol
754    end if
755
756    ! Locate the surface half-level
757    if (pressure_hl(1,1) > pressure_hl(1,2)) then
758      isurf = 1
759    else
760      isurf = nlev+1
761    end if
762
763    do jcol = i1,i2
764      eta = (pressure_hl(jcol,1:nlev)+pressure_hl(jcol,2:nlev+1)) &
765           &  * (0.5_jprb / pressure_hl(jcol,isurf))
766      eff_separation = coeff_a + coeff_b * exp(-eta**power)
767      this%inv_cloud_effective_size(jcol,:) = 1.0_jprb / (eff_separation &
768           &  * sqrt(max(1.0e-5_jprb,this%fraction(jcol,:)*(1.0_jprb-this%fraction(jcol,:)))))
769      this%inv_inhom_effective_size(jcol,:) = 1.0_jprb / (eff_separation * inhom_sep_factor &
770           &  * sqrt(max(1.0e-5_jprb,0.5_jprb*this%fraction(jcol,:)*(1.0_jprb-0.5_jprb*this%fraction(jcol,:)))))
771    end do
772
773    if (lhook) call dr_hook('radiation_cloud:param_cloud_effective_separation_eta',1,hook_handle)
774
775  end subroutine param_cloud_effective_separation_eta
776
777
778  !---------------------------------------------------------------------
779  ! Remove "ghost" clouds: those with a cloud fraction that is too
780  ! small to treat sensibly (e.g. because it implies that the
781  ! "in-cloud" water content is too high), or with a cloud water
782  ! content that is too small.  We do this in one place to ensure that
783  ! all subsequent subroutines can assume that if cloud_fraction > 0.0
784  ! then cloud is really present and should be treated.
785  subroutine crop_cloud_fraction(this, istartcol, iendcol, &
786       &    cloud_fraction_threshold, cloud_mixing_ratio_threshold)
787   
788    use yomhook, only : lhook, dr_hook, jphook
789
790    class(cloud_type), intent(inout) :: this
791    integer,           intent(in)    :: istartcol, iendcol
792
793    integer :: nlev, ntype
794    integer :: jcol, jlev, jh
795
796    real(jprb) :: cloud_fraction_threshold, cloud_mixing_ratio_threshold
797    real(jprb) :: sum_mixing_ratio(istartcol:iendcol)
798
799    real(jphook) :: hook_handle
800
801    if (lhook) call dr_hook('radiation_cloud:crop_cloud_fraction',0,hook_handle)
802
803    nlev  = size(this%fraction,2)
804    ntype = size(this%mixing_ratio,3)
805   
806    do jlev = 1,nlev
807      do jcol = istartcol,iendcol
808        sum_mixing_ratio(jcol) = 0.0_jprb
809      end do
810      do jh = 1, ntype
811        do jcol = istartcol,iendcol
812          sum_mixing_ratio(jcol) = sum_mixing_ratio(jcol) + this%mixing_ratio(jcol,jlev,jh)
813        end do
814      end do
815      do jcol = istartcol,iendcol
816        if (this%fraction(jcol,jlev)        < cloud_fraction_threshold &
817             &  .or. sum_mixing_ratio(jcol) < cloud_mixing_ratio_threshold) then
818          this%fraction(jcol,jlev) = 0.0_jprb
819        end if
820      end do
821    end do
822
823    if (lhook) call dr_hook('radiation_cloud:crop_cloud_fraction',1,hook_handle)
824
825  end subroutine crop_cloud_fraction
826
827
828  !---------------------------------------------------------------------
829  ! Return .true. if variables are out of a physically sensible range,
830  ! optionally only considering columns between istartcol and iendcol
831  function out_of_physical_bounds(this, istartcol, iendcol, do_fix) result(is_bad)
832
833    use yomhook,          only : lhook, dr_hook, jphook
834    use radiation_check, only : out_of_bounds_2d, out_of_bounds_3d
835
836    class(cloud_type), intent(inout) :: this
837    integer,  optional,intent(in) :: istartcol, iendcol
838    logical,  optional,intent(in) :: do_fix
839    logical                       :: is_bad
840
841    logical    :: do_fix_local
842
843    real(jphook) :: hook_handle
844
845    if (lhook) call dr_hook('radiation_cloud:out_of_physical_bounds',0,hook_handle)
846
847    if (present(do_fix)) then
848      do_fix_local = do_fix
849    else
850      do_fix_local = .false.
851    end if
852
853    is_bad =    out_of_bounds_3d(this%mixing_ratio, 'cloud%mixing_ratio', 0.0_jprb, 1.0_jprb, &
854         &                       do_fix_local, i1=istartcol, i2=iendcol) &
855         & .or. out_of_bounds_3d(this%effective_radius, 'cloud%effective_radius', 0.0_jprb, 0.1_jprb, &
856         &                       do_fix_local, i1=istartcol, i2=iendcol) &
857         & .or. out_of_bounds_2d(this%fraction, 'cloud%fraction', 0.0_jprb, 1.0_jprb, &
858         &                       do_fix_local, i1=istartcol, i2=iendcol) &
859         & .or. out_of_bounds_2d(this%fractional_std, 'fractional_std', 0.0_jprb, 10.0_jprb, &
860         &                       do_fix_local, i1=istartcol, i2=iendcol) &
861         & .or. out_of_bounds_2d(this%inv_cloud_effective_size, 'inv_cloud_effective_size', &
862         &                       0.0_jprb, 1.0_jprb, do_fix_local, i1=istartcol, i2=iendcol) &
863         & .or. out_of_bounds_2d(this%inv_inhom_effective_size, 'inv_inhom_effective_size', &
864         &                       0.0_jprb, 1.0_jprb, do_fix_local, i1=istartcol, i2=iendcol) &
865         & .or. out_of_bounds_2d(this%overlap_param, 'overlap_param', -0.5_jprb, 1.0_jprb, &
866         &                       do_fix_local, i1=istartcol, i2=iendcol)
867
868    if (lhook) call dr_hook('radiation_cloud:out_of_physical_bounds',1,hook_handle)
869
870  end function out_of_physical_bounds
871 
872end module radiation_cloud
Note: See TracBrowser for help on using the repository browser.