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

Last change on this file since 5690 was 5690, checked in by lguez, 39 hours ago

Do not allocate inv_cloud_effective_size

This field was not defined if Spartacus was not used so it appeared
with random values in file inputs.nc created by procedure
save_inputs. It is allocated elsewhere if Spartacus is used. We are
cherry-picking commit c7a36b61 in ecrad repository. See pull requests
35 and 18 in upstream ECRad.

File size: 34.6 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
145    if (present(use_inhom_effective_size)) then
146      if (use_inhom_effective_size) then
147        allocate(this%inv_inhom_effective_size(ncol,nlev))
148      end if
149    end if
150
151    if (lhook) call dr_hook('radiation_cloud:allocate',1,hook_handle)
152
153  end subroutine allocate_cloud_arrays
154
155
156  !---------------------------------------------------------------------
157  ! Deallocate arrays
158  subroutine deallocate_cloud_arrays(this)
159
160    use yomhook,     only : lhook, dr_hook, jphook
161
162    class(cloud_type), intent(inout) :: this
163
164    real(jphook) :: hook_handle
165
166    if (lhook) call dr_hook('radiation_cloud:deallocate',0,hook_handle)
167
168    nullify(this%q_liq)
169    nullify(this%q_ice)
170    nullify(this%re_liq)
171    nullify(this%re_ice)
172
173    if (allocated(this%mixing_ratio))     deallocate(this%mixing_ratio)
174    if (allocated(this%effective_radius)) deallocate(this%effective_radius)
175    if (allocated(this%fraction))         deallocate(this%fraction)
176    if (allocated(this%overlap_param))    deallocate(this%overlap_param)
177    if (allocated(this%fractional_std))   deallocate(this%fractional_std)
178    if (allocated(this%inv_cloud_effective_size)) &
179         &  deallocate(this%inv_cloud_effective_size)
180    if (allocated(this%inv_inhom_effective_size)) &
181         &  deallocate(this%inv_inhom_effective_size)
182
183    if (lhook) call dr_hook('radiation_cloud:deallocate',1,hook_handle)
184
185  end subroutine deallocate_cloud_arrays
186
187
188  !---------------------------------------------------------------------
189  ! Compute and store the overlap parameter from the provided overlap
190  ! decorrelation length (in metres).  If istartcol and/or iendcol are
191  ! provided then only columns in this range are computed.  If the
192  ! overlap_param array has not been allocated then it will be
193  ! allocated to be of the correct size relative to the pressure
194  ! field. This version assumes a fixed decorrelation_length for all
195  ! columns.
196  subroutine set_overlap_param_fix(this, thermodynamics, decorrelation_length, &
197       &  istartcol, iendcol)
198
199    use yomhook,                  only : lhook, dr_hook, jphook
200    use radiation_thermodynamics, only : thermodynamics_type
201    use radiation_constants,      only : GasConstantDryAir, AccelDueToGravity
202
203    class(cloud_type),         intent(inout) :: this
204    type(thermodynamics_type), intent(in)    :: thermodynamics
205    real(jprb),                intent(in)    :: decorrelation_length ! m
206    integer,         optional, intent(in)    :: istartcol, iendcol
207
208    ! Ratio of gas constant for dry air to acceleration due to gravity
209    real(jprb), parameter :: R_over_g = GasConstantDryAir / AccelDueToGravity
210
211    ! Process only columns i1 to i2, which will be istartcol to
212    ! iendcol if they were provided
213    integer :: i1, i2
214
215    integer :: ncol, nlev
216
217    integer :: jcol, jlev
218
219    real(jphook) :: hook_handle
220
221    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_fix',0,hook_handle)
222
223    ! Pressure at half-levels, pressure_hl, is defined at nlev+1
224    ! points
225    ncol = size(thermodynamics%pressure_hl,dim=1)
226    nlev = size(thermodynamics%pressure_hl,dim=2)-1
227
228    if (present(istartcol)) then
229      i1 = istartcol
230    else
231      i1 = 1
232    end if
233
234    if (present(iendcol)) then
235      i2 = iendcol
236    else
237      i2 = ncol
238    end if
239
240    if (.not. allocated(this%overlap_param)) then
241      ! If pressure is of size (ncol,nlev+1) then overlap_param is of
242      ! size (ncol,nlev-1), since overlap parameter is only defined here
243      ! for interfaces between model layers, not for the interface to
244      ! space or the surface
245      allocate(this%overlap_param(ncol, nlev-1))
246    end if
247
248    if (thermodynamics%pressure_hl(i1,2) > thermodynamics%pressure_hl(i1,1)) then
249      ! Pressure is increasing with index (order of layers is
250      ! top-of-atmosphere to surface). In case pressure_hl(:,1)=0, we
251      ! don't take the logarithm of the first pressure in each column.
252      do jcol = i1,i2
253        this%overlap_param(jcol,1) = exp(-(R_over_g/decorrelation_length) &
254             &                            * thermodynamics%temperature_hl(jcol,2) &
255             &                            *log(thermodynamics%pressure_hl(jcol,3) &
256             &                                /thermodynamics%pressure_hl(jcol,2)))
257      end do
258
259      do jlev = 2,nlev-1
260        do jcol = i1,i2
261          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length) &
262              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
263              &                            *log(thermodynamics%pressure_hl(jcol,jlev+2) &
264              &                                /thermodynamics%pressure_hl(jcol,jlev)))
265        end do
266      end do
267
268    else
269       ! Pressure is decreasing with index (order of layers is surface
270       ! to top-of-atmosphere).  In case pressure_hl(:,nlev+1)=0, we
271       ! don't take the logarithm of the last pressure in each column.
272      do jlev = 1,nlev-2
273        do jcol = i1,i2
274          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length) &
275              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
276              &                            *log(thermodynamics%pressure_hl(jcol,jlev) &
277              &                                /thermodynamics%pressure_hl(jcol,jlev+2)))
278        end do
279      end do
280
281      do jcol = i1,i2
282        this%overlap_param(jcol,nlev-1) = exp(-(R_over_g/decorrelation_length) &
283            &                            * thermodynamics%temperature_hl(jcol,nlev) &
284            &                            *log(thermodynamics%pressure_hl(jcol,nlev-1) &
285            &                                /thermodynamics%pressure_hl(jcol,nlev)))
286      end do
287    end if
288
289    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_fix',1,hook_handle)
290
291  end subroutine set_overlap_param_fix
292
293
294  !---------------------------------------------------------------------
295  ! Compute and store the overlap parameter from the provided overlap
296  ! decorrelation length (in metres), which may vary with column. Only
297  ! columns from istartcol to iendcol are computed.  If the
298  ! overlap_param array has not been allocated then it will be
299  ! allocated to be of the correct size relative to the pressure
300  ! field.
301  subroutine set_overlap_param_var(this, thermodynamics, decorrelation_length, &
302       &                           istartcol, iendcol)
303
304    use yomhook,                  only : lhook, dr_hook, jphook
305    use radiation_thermodynamics, only : thermodynamics_type
306    use radiation_constants,      only : GasConstantDryAir, AccelDueToGravity
307
308    class(cloud_type),         intent(inout) :: this
309    type(thermodynamics_type), intent(in)    :: thermodynamics
310    integer,                   intent(in)    :: istartcol, iendcol
311    real(jprb),                intent(in)    :: decorrelation_length(istartcol:iendcol) ! m
312
313    ! Ratio of gas constant for dry air to acceleration due to gravity
314    real(jprb), parameter :: R_over_g = GasConstantDryAir / AccelDueToGravity
315
316    integer :: ncol, nlev
317
318    integer :: jcol, jlev
319
320    real(jphook) :: hook_handle
321
322    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_var',0,hook_handle)
323
324    ! Pressure at half-levels, pressure_hl, is defined at nlev+1
325    ! points
326    ncol = size(thermodynamics%pressure_hl,dim=1)
327    nlev = size(thermodynamics%pressure_hl,dim=2)-1
328
329    if (.not. allocated(this%overlap_param)) then
330      ! If pressure is of size (ncol,nlev+1) then overlap_param is of
331      ! size (ncol,nlev-1), since overlap parameter is only defined here
332      ! for interfaces between model layers, not for the interface to
333      ! space or the surface
334      allocate(this%overlap_param(ncol, nlev-1))
335    end if
336
337    if (thermodynamics%pressure_hl(istartcol,2) > thermodynamics%pressure_hl(istartcol,1)) then
338      ! Pressure is increasing with index (order of layers is
339      ! top-of-atmosphere to surface). In case pressure_hl(:,1)=0, we
340      ! don't take the logarithm of the first pressure in each column.
341      do jcol = istartcol,iendcol
342        this%overlap_param(jcol,1) = exp(-(R_over_g/decorrelation_length(jcol)) &
343             &                            * thermodynamics%temperature_hl(jcol,2) &
344             &                            *log(thermodynamics%pressure_hl(jcol,3) &
345             &                                /thermodynamics%pressure_hl(jcol,2)))
346      end do
347
348      do jlev = 2,nlev-1
349        do jcol = istartcol,iendcol
350          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length(jcol)) &
351              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
352              &                            *log(thermodynamics%pressure_hl(jcol,jlev+2) &
353              &                                /thermodynamics%pressure_hl(jcol,jlev)))
354        end do
355      end do
356
357    else
358       ! Pressure is decreasing with index (order of layers is surface
359       ! to top-of-atmosphere).  In case pressure_hl(:,nlev+1)=0, we
360       ! don't take the logarithm of the last pressure in each column.
361      do jlev = 1,nlev-2
362        do jcol = istartcol,iendcol
363          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length(jcol)) &
364              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
365              &                            *log(thermodynamics%pressure_hl(jcol,jlev) &
366              &                                /thermodynamics%pressure_hl(jcol,jlev+2)))
367        end do
368      end do
369
370      do jcol = istartcol,iendcol
371        this%overlap_param(jcol,nlev-1) = exp(-(R_over_g/decorrelation_length(jcol)) &
372            &                            * thermodynamics%temperature_hl(jcol,nlev) &
373            &                            *log(thermodynamics%pressure_hl(jcol,nlev-1) &
374            &                                /thermodynamics%pressure_hl(jcol,nlev)))
375      end do
376    end if
377
378    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_var',1,hook_handle)
379
380  end subroutine set_overlap_param_var
381
382! AI 04 2024 variation de la longueur Ld en fonction de la verticale 
383  subroutine set_overlap_param_var2D(this, thermodynamics, decorrelation_length, &
384       &                           klev, istartcol, iendcol)
385
386    use yomhook,                  only : lhook, dr_hook, jphook
387    use radiation_thermodynamics, only : thermodynamics_type
388    use radiation_constants,      only : GasConstantDryAir, AccelDueToGravity
389
390    integer,                   intent(in)    :: klev
391    class(cloud_type),         intent(inout) :: this
392    type(thermodynamics_type), intent(in)    :: thermodynamics
393    integer,                   intent(in)    :: istartcol, iendcol
394    real(jprb),                intent(in)    :: decorrelation_length(istartcol:iendcol,klev) ! m
395
396    ! Ratio of gas constant for dry air to acceleration due to gravity
397    real(jprb), parameter :: R_over_g = GasConstantDryAir / AccelDueToGravity
398
399    integer :: ncol, nlev
400
401    integer :: jcol, jlev
402
403    real(jphook) :: hook_handle
404
405    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_var',0,hook_handle)
406
407    ! Pressure at half-levels, pressure_hl, is defined at nlev+1
408    ! points
409    ncol = size(thermodynamics%pressure_hl,dim=1)
410    nlev = size(thermodynamics%pressure_hl,dim=2)-1
411
412    if (.not. allocated(this%overlap_param)) then
413      ! If pressure is of size (ncol,nlev+1) then overlap_param is of
414      ! size (ncol,nlev-1), since overlap parameter is only defined here
415      ! for interfaces between model layers, not for the interface to
416      ! space or the surface
417      allocate(this%overlap_param(ncol, nlev-1))
418    end if
419
420    if (thermodynamics%pressure_hl(istartcol,2) > thermodynamics%pressure_hl(istartcol,1)) then
421      ! Pressure is increasing with index (order of layers is
422      ! top-of-atmosphere to surface). In case pressure_hl(:,1)=0, we
423      ! don't take the logarithm of the first pressure in each column.
424      do jcol = istartcol,iendcol
425        this%overlap_param(jcol,1) = exp(-(R_over_g/decorrelation_length(jcol,1)) &
426             &                            * thermodynamics%temperature_hl(jcol,2) &
427             &                            *log(thermodynamics%pressure_hl(jcol,3) &
428             &                                /thermodynamics%pressure_hl(jcol,2)))
429      end do
430
431      do jlev = 2,nlev-1
432        do jcol = istartcol,iendcol
433          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length(jcol,jlev)) &
434              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
435              &                            *log(thermodynamics%pressure_hl(jcol,jlev+2) &
436              &                                /thermodynamics%pressure_hl(jcol,jlev)))
437        end do
438      end do
439
440    else
441       ! Pressure is decreasing with index (order of layers is surface
442       ! to top-of-atmosphere).  In case pressure_hl(:,nlev+1)=0, we
443       ! don't take the logarithm of the last pressure in each column.
444      do jlev = 1,nlev-2
445        do jcol = istartcol,iendcol
446          this%overlap_param(jcol,jlev) = exp(-(0.5_jprb*R_over_g/decorrelation_length(jcol,jlev)) &
447              &                            * thermodynamics%temperature_hl(jcol,jlev+1) &
448              &                            *log(thermodynamics%pressure_hl(jcol,jlev) &
449              &                                /thermodynamics%pressure_hl(jcol,jlev+2)))
450        end do
451      end do
452
453      do jcol = istartcol,iendcol
454      ! AI ATTENTION a verifier decorrelation_length(jcol,nlev-1) ou nlev
455        this%overlap_param(jcol,nlev-1) = exp(-(R_over_g/decorrelation_length(jcol,nlev-1)) &
456            &                            * thermodynamics%temperature_hl(jcol,nlev) &
457            &                            *log(thermodynamics%pressure_hl(jcol,nlev-1) &
458            &                                /thermodynamics%pressure_hl(jcol,nlev)))
459      end do
460    end if
461
462    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_var',1,hook_handle)
463
464  end subroutine set_overlap_param_var2D
465
466  !---------------------------------------------------------------------
467  ! Compute and store the overlap parameter from the provided overlap
468  ! decorrelation length (in metres).  If istartcol and/or iendcol are
469  ! provided then only columns in this range are computed.  If the
470  ! overlap_param array has not been allocated then it will be
471  ! allocated to be of the correct size relative to the pressure
472  ! field. This is the APPROXIMATE method as it assumes a fixed
473  ! atmospheric scale height, which leads to differences particularly
474  ! in low cloud.
475  subroutine set_overlap_param_approx(this, thermodynamics, decorrelation_length, &
476       &  istartcol, iendcol)
477
478    use yomhook,                  only : lhook, dr_hook, jphook
479    use radiation_thermodynamics, only : thermodynamics_type
480
481    class(cloud_type),         intent(inout) :: this
482    type(thermodynamics_type), intent(in)    :: thermodynamics
483    real(jprb),                intent(in)    :: decorrelation_length ! m
484    integer,         optional, intent(in)    :: istartcol, iendcol
485
486    ! To convert decorrelation length (m) to overlap parameter between
487    ! layers, we need an estimate for the thickness of the layer. This
488    ! is found using the pressure difference between the edges of the
489    ! layer, along with the approximate scale height of the atmosphere
490    ! (m) given here:
491    real(jprb), parameter :: scale_height = 8000.0_jprb
492
493    ! Process only columns i1 to i2, which will be istartcol to
494    ! iendcol if they were provided
495    integer :: i1, i2
496
497    integer :: ncol, nlev
498
499    real(jphook) :: hook_handle
500
501    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_approx',0,hook_handle)
502
503    ! Pressure at half-levels, pressure_hl, is defined at nlev+1
504    ! points
505    ncol = size(thermodynamics%pressure_hl,dim=1)
506    nlev = size(thermodynamics%pressure_hl,dim=2)-1
507
508    if (present(istartcol)) then
509      i1 = istartcol
510    else
511      i1 = 1
512    end if
513
514    if (present(iendcol)) then
515      i2 = iendcol
516    else
517      i2 = ncol
518    end if
519
520    if (.not. allocated(this%overlap_param)) then
521      ! If pressure is of size (ncol,nlev+1) then overlap_param is of
522      ! size (ncol,nlev-1), since overlap parameter is only defined here
523      ! for interfaces between model layers, not for the interface to
524      ! space or the surface
525      allocate(this%overlap_param(ncol, nlev-1))
526    end if
527
528    if (thermodynamics%pressure_hl(i1,2) > thermodynamics%pressure_hl(i1,1)) then
529       ! Pressure is increasing with index (order of layers is
530       ! top-of-atmosphere to surface). In case pressure_hl(:,1)=0, we
531       ! don't take the logarithm of the first pressure in each
532       ! column.
533       this%overlap_param(i1:i2,:) = exp(-(scale_height/decorrelation_length) &
534            &  * ( log(thermodynamics%pressure_hl(i1:i2,3:nlev+1) &
535            &         /thermodynamics%pressure_hl(i1:i2,2:nlev  )) ) )
536    else
537       ! Pressure is decreasing with index (order of layers is surface
538       ! to top-of-atmosphere).  In case pressure_hl(:,nlev+1)=0, we
539       ! don't take the logarithm of the last pressure in each column.
540       this%overlap_param(i1:i2,:) = exp(-(scale_height/decorrelation_length) &
541            &  * ( log(thermodynamics%pressure_hl(i1:i2,1:nlev-1) &
542            &         /thermodynamics%pressure_hl(i1:i2,2:nlev  )) ) )
543    end if
544
545    if (lhook) call dr_hook('radiation_cloud:set_overlap_param_approx',1,hook_handle)
546
547  end subroutine set_overlap_param_approx
548
549
550  !---------------------------------------------------------------------
551  ! Create a matrix of constant fractional standard deviations
552  ! (dimensionless)
553  subroutine create_fractional_std(this, ncol, nlev, frac_std)
554
555    use yomhook,                  only : lhook, dr_hook, jphook
556
557    class(cloud_type), intent(inout) :: this
558    integer,           intent(in)    :: ncol, nlev
559    real(jprb),        intent(in)    :: frac_std
560
561    real(jphook) :: hook_handle
562
563    if (lhook) call dr_hook('radiation_cloud:create_fractional_std',0,hook_handle)
564
565    if (allocated(this%fractional_std)) then
566       deallocate(this%fractional_std)
567    end if
568   
569    allocate(this%fractional_std(ncol, nlev))
570
571    this%fractional_std = frac_std
572
573    if (lhook) call dr_hook('radiation_cloud:create_fractional_std',1,hook_handle)
574
575  end subroutine create_fractional_std
576
577
578  !---------------------------------------------------------------------
579  ! Create a matrix of constant inverse cloud effective size (m-1)
580  subroutine create_inv_cloud_effective_size(this, ncol, nlev, inv_eff_size)
581
582    use yomhook,                  only : lhook, dr_hook, jphook
583
584    class(cloud_type), intent(inout) :: this
585    integer,           intent(in)    :: ncol, nlev
586    real(jprb),        intent(in)    :: inv_eff_size
587
588    real(jphook) :: hook_handle
589
590    if (lhook) call dr_hook('radiation_cloud:create_inv_cloud_effective_size',0,hook_handle)
591
592    if (allocated(this%inv_cloud_effective_size)) then
593       deallocate(this%inv_cloud_effective_size)
594    end if
595   
596    allocate(this%inv_cloud_effective_size(ncol, nlev))
597
598    this%inv_cloud_effective_size = inv_eff_size
599
600    if (lhook) call dr_hook('radiation_cloud:create_inv_cloud_effective_size',1,hook_handle)
601
602  end subroutine create_inv_cloud_effective_size
603
604
605  !---------------------------------------------------------------------
606  ! Create a matrix of inverse cloud effective size (m-1) according to
607  ! the value of eta (=pressure divided by surface pressure)
608  subroutine create_inv_cloud_effective_size_eta(this, ncol, nlev, &
609       &  pressure_hl, inv_eff_size_low, inv_eff_size_mid, inv_eff_size_high, &
610       &  eta_low_mid, eta_mid_high, istartcol, iendcol)
611
612    use yomhook,                  only : lhook, dr_hook, jphook
613
614    class(cloud_type), intent(inout) :: this
615    integer,           intent(in)    :: ncol, nlev
616    ! Pressure on half levels (Pa)
617    real(jprb),        intent(in)    :: pressure_hl(:,:)
618    ! Inverse effective size for low, mid and high cloud (m-1)
619    real(jprb),        intent(in)    :: inv_eff_size_low
620    real(jprb),        intent(in)    :: inv_eff_size_mid
621    real(jprb),        intent(in)    :: inv_eff_size_high
622    ! Eta values at low-mid and mid-high interfaces
623    real(jprb),        intent(in)    :: eta_low_mid, eta_mid_high
624    integer, optional, intent(in)    :: istartcol, iendcol
625
626    ! Ratio of layer midpoint pressure to surface pressure
627    real(jprb) :: eta(nlev)
628
629    ! Indices of column, level and surface half-level
630    integer :: jcol, isurf
631
632    ! Local values of istartcol, iendcol
633    integer :: i1, i2
634
635    real(jphook) :: hook_handle
636
637    if (lhook) call dr_hook('radiation_cloud:create_inv_cloud_effective_size_eta',0,hook_handle)
638
639    if (allocated(this%inv_cloud_effective_size)) then
640      deallocate(this%inv_cloud_effective_size)
641    end if
642   
643    allocate(this%inv_cloud_effective_size(ncol, nlev))
644
645    if (present(istartcol)) then
646      i1 = istartcol
647    else
648      i1 = 1
649    end if
650
651    if (present(iendcol)) then
652      i2 = iendcol
653    else
654      i2 = ncol
655    end if
656
657    ! Locate the surface half-level
658    if (pressure_hl(1,1) > pressure_hl(1,2)) then
659      isurf = 1
660    else
661      isurf = nlev+1
662    end if
663
664    do jcol = i1,i2
665      eta = (pressure_hl(jcol,1:nlev)+pressure_hl(jcol,2:nlev+1)) &
666           &  * (0.5_jprb / pressure_hl(jcol,isurf))
667      where (eta > eta_low_mid)
668        this%inv_cloud_effective_size(jcol,:) = inv_eff_size_low
669      elsewhere (eta > eta_mid_high)
670        this%inv_cloud_effective_size(jcol,:) = inv_eff_size_mid
671      elsewhere
672        this%inv_cloud_effective_size(jcol,:) = inv_eff_size_high
673      end where
674    end do
675
676    if (lhook) call dr_hook('radiation_cloud:create_inv_cloud_effective_size_eta',1,hook_handle)
677
678  end subroutine create_inv_cloud_effective_size_eta
679
680
681  !---------------------------------------------------------------------
682  ! Create a matrix of inverse cloud and inhomogeneity effective size
683  ! (m-1) parameterized according to the value of eta (=pressure
684  ! divided by surface pressure): effective_separation =
685  ! coeff_a + coeff_b*exp(-(eta**power)). 
686  subroutine param_cloud_effective_separation_eta(this, ncol, nlev, &
687       &  pressure_hl, separation_surf, separation_toa, power, &
688       &  inhom_separation_factor, istartcol, iendcol)
689
690    use yomhook,                  only : lhook, dr_hook, jphook
691
692    class(cloud_type), intent(inout) :: this
693    integer,           intent(in)    :: ncol, nlev
694    ! Pressure on half levels (Pa)
695    real(jprb),        intent(in)    :: pressure_hl(:,:)
696    ! Separation distances at surface and top-of-atmosphere, and power
697    ! on eta
698    real(jprb),           intent(in) :: separation_surf ! m
699    real(jprb),           intent(in) :: separation_toa ! m
700    real(jprb),           intent(in) :: power
701    real(jprb), optional, intent(in) :: inhom_separation_factor
702    integer,    optional, intent(in) :: istartcol, iendcol
703
704    ! Ratio of layer midpoint pressure to surface pressure
705    real(jprb) :: eta(nlev)
706
707    ! Effective cloud separation (m)
708    real(jprb) :: eff_separation(nlev)
709
710    ! Coefficients used to compute effective separation distance
711    real(jprb) :: coeff_e, coeff_a, coeff_b, inhom_sep_factor
712
713    ! Indices of column, level and surface half-level
714    integer :: jcol, isurf
715
716    ! Local values of istartcol, iendcol
717    integer :: i1, i2
718
719    real(jphook) :: hook_handle
720
721    if (lhook) call dr_hook('radiation_cloud:param_cloud_effective_separation_eta',0,hook_handle)
722
723    if (present(inhom_separation_factor)) then
724      inhom_sep_factor = inhom_separation_factor
725    else
726      inhom_sep_factor = 1.0_jprb
727    end if
728
729    coeff_e = 1.0_jprb - exp(-1.0_jprb)
730    coeff_b = (separation_toa - separation_surf) / coeff_e
731    coeff_a = separation_toa - coeff_b
732
733    if (allocated(this%inv_cloud_effective_size)) then
734      deallocate(this%inv_cloud_effective_size)
735    end if
736     if (allocated(this%inv_inhom_effective_size)) then
737      deallocate(this%inv_inhom_effective_size)
738    end if
739   
740    allocate(this%inv_cloud_effective_size(ncol, nlev))
741    allocate(this%inv_inhom_effective_size(ncol, nlev))
742
743    if (present(istartcol)) then
744      i1 = istartcol
745    else
746      i1 = 1
747    end if
748
749    if (present(iendcol)) then
750      i2 = iendcol
751    else
752      i2 = ncol
753    end if
754
755    ! Locate the surface half-level
756    if (pressure_hl(1,1) > pressure_hl(1,2)) then
757      isurf = 1
758    else
759      isurf = nlev+1
760    end if
761
762    do jcol = i1,i2
763      eta = (pressure_hl(jcol,1:nlev)+pressure_hl(jcol,2:nlev+1)) &
764           &  * (0.5_jprb / pressure_hl(jcol,isurf))
765      eff_separation = coeff_a + coeff_b * exp(-eta**power)
766      this%inv_cloud_effective_size(jcol,:) = 1.0_jprb / (eff_separation &
767           &  * sqrt(max(1.0e-5_jprb,this%fraction(jcol,:)*(1.0_jprb-this%fraction(jcol,:)))))
768      this%inv_inhom_effective_size(jcol,:) = 1.0_jprb / (eff_separation * inhom_sep_factor &
769           &  * sqrt(max(1.0e-5_jprb,0.5_jprb*this%fraction(jcol,:)*(1.0_jprb-0.5_jprb*this%fraction(jcol,:)))))
770    end do
771
772    if (lhook) call dr_hook('radiation_cloud:param_cloud_effective_separation_eta',1,hook_handle)
773
774  end subroutine param_cloud_effective_separation_eta
775
776
777  !---------------------------------------------------------------------
778  ! Remove "ghost" clouds: those with a cloud fraction that is too
779  ! small to treat sensibly (e.g. because it implies that the
780  ! "in-cloud" water content is too high), or with a cloud water
781  ! content that is too small.  We do this in one place to ensure that
782  ! all subsequent subroutines can assume that if cloud_fraction > 0.0
783  ! then cloud is really present and should be treated.
784  subroutine crop_cloud_fraction(this, istartcol, iendcol, &
785       &    cloud_fraction_threshold, cloud_mixing_ratio_threshold)
786   
787    use yomhook, only : lhook, dr_hook, jphook
788
789    class(cloud_type), intent(inout) :: this
790    integer,           intent(in)    :: istartcol, iendcol
791
792    integer :: nlev, ntype
793    integer :: jcol, jlev, jh
794
795    real(jprb) :: cloud_fraction_threshold, cloud_mixing_ratio_threshold
796    real(jprb) :: sum_mixing_ratio(istartcol:iendcol)
797
798    real(jphook) :: hook_handle
799
800    if (lhook) call dr_hook('radiation_cloud:crop_cloud_fraction',0,hook_handle)
801
802    nlev  = size(this%fraction,2)
803    ntype = size(this%mixing_ratio,3)
804   
805    do jlev = 1,nlev
806      do jcol = istartcol,iendcol
807        sum_mixing_ratio(jcol) = 0.0_jprb
808      end do
809      do jh = 1, ntype
810        do jcol = istartcol,iendcol
811          sum_mixing_ratio(jcol) = sum_mixing_ratio(jcol) + this%mixing_ratio(jcol,jlev,jh)
812        end do
813      end do
814      do jcol = istartcol,iendcol
815        if (this%fraction(jcol,jlev)        < cloud_fraction_threshold &
816             &  .or. sum_mixing_ratio(jcol) < cloud_mixing_ratio_threshold) then
817          this%fraction(jcol,jlev) = 0.0_jprb
818        end if
819      end do
820    end do
821
822    if (lhook) call dr_hook('radiation_cloud:crop_cloud_fraction',1,hook_handle)
823
824  end subroutine crop_cloud_fraction
825
826
827  !---------------------------------------------------------------------
828  ! Return .true. if variables are out of a physically sensible range,
829  ! optionally only considering columns between istartcol and iendcol
830  function out_of_physical_bounds(this, istartcol, iendcol, do_fix) result(is_bad)
831
832    use yomhook,          only : lhook, dr_hook, jphook
833    use radiation_check, only : out_of_bounds_2d, out_of_bounds_3d
834
835    class(cloud_type), intent(inout) :: this
836    integer,  optional,intent(in) :: istartcol, iendcol
837    logical,  optional,intent(in) :: do_fix
838    logical                       :: is_bad
839
840    logical    :: do_fix_local
841
842    real(jphook) :: hook_handle
843
844    if (lhook) call dr_hook('radiation_cloud:out_of_physical_bounds',0,hook_handle)
845
846    if (present(do_fix)) then
847      do_fix_local = do_fix
848    else
849      do_fix_local = .false.
850    end if
851
852    is_bad =    out_of_bounds_3d(this%mixing_ratio, 'cloud%mixing_ratio', 0.0_jprb, 1.0_jprb, &
853         &                       do_fix_local, i1=istartcol, i2=iendcol) &
854         & .or. out_of_bounds_3d(this%effective_radius, 'cloud%effective_radius', 0.0_jprb, 0.1_jprb, &
855         &                       do_fix_local, i1=istartcol, i2=iendcol) &
856         & .or. out_of_bounds_2d(this%fraction, 'cloud%fraction', 0.0_jprb, 1.0_jprb, &
857         &                       do_fix_local, i1=istartcol, i2=iendcol) &
858         & .or. out_of_bounds_2d(this%fractional_std, 'fractional_std', 0.0_jprb, 10.0_jprb, &
859         &                       do_fix_local, i1=istartcol, i2=iendcol) &
860         & .or. out_of_bounds_2d(this%inv_cloud_effective_size, 'inv_cloud_effective_size', &
861         &                       0.0_jprb, 1.0_jprb, do_fix_local, i1=istartcol, i2=iendcol) &
862         & .or. out_of_bounds_2d(this%inv_inhom_effective_size, 'inv_inhom_effective_size', &
863         &                       0.0_jprb, 1.0_jprb, do_fix_local, i1=istartcol, i2=iendcol) &
864         & .or. out_of_bounds_2d(this%overlap_param, 'overlap_param', -0.5_jprb, 1.0_jprb, &
865         &                       do_fix_local, i1=istartcol, i2=iendcol)
866
867    if (lhook) call dr_hook('radiation_cloud:out_of_physical_bounds',1,hook_handle)
868
869  end function out_of_physical_bounds
870 
871end module radiation_cloud
Note: See TracBrowser for help on using the repository browser.