source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_cloud_cover.F90 @ 4911

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

Ecrad update in LMDZ / Implementation of Ecrad double call in LMDZ

  • version 1.6.1 (commit 210d7911380f53a788c3cad73b3cf9b4e022ef87)
  • interface routines between lmdz and ecrad grouped in a new "lmdz" directory
  • double call method redesigned so as to go through the Ecrad initialization and configuration part only once for the entire simulation
  • clean-up in the read.F routine: delete unitules arguments
  • modification of default gas model in namelist (default: ECCKD)
File size: 19.5 KB
Line 
1! radiation_cloud_cover.F90 - Compute cumulative cloud cover for McICA
2!
3! (C) Copyright 2016- 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! Generate profiles of the cumulative cloud cover as seen from TOA,
16! used in the McICA cloud generator.
17!
18! Modifications
19!   2020-10-07  R. Hogan  Ensure iobj1 initialized in case of alpha_obj==0
20
21#include "ecrad_config.h"
22
23module radiation_cloud_cover
24
25  use parkind1, only           : jprb
26
27  public
28
29  ! Three overlap schemes.  Note that "Exponential" means that
30  ! clear-sky regions have no special significance for computing the
31  ! cumulative cloud cover: non-contiguous clouds are exponentially
32  ! rather than randomly overlapped. This is the situaition in the
33  ! McRad radiation scheme at ECMWF.
34  enum, bind(c)
35    enumerator IOverlapMaximumRandom, IOverlapExponentialRandom, &
36         &     IOverlapExponential
37  end enum
38  character(len=*), parameter :: OverlapName(0:2) = (/ 'Max-Ran', &
39       &                                               'Exp-Ran', &
40       &                                               'Exp-Exp' /)
41
42  ! Maximum cloud fraction distinguishable from 1
43  real(jprb), parameter :: MaxCloudFrac = 1.0_jprb-epsilon(1.0_jprb)*10.0_jprb
44
45
46contains
47
48  !---------------------------------------------------------------------
49  ! Convert "beta" overlap parameter of Shonk et al. (2010) to "alpha"
50  ! overlap parameter of Hogan and Illingworth (2000)
51  elemental function beta2alpha(beta, frac1, frac2)
52
53    implicit none
54
55    ! Beta overlap parameter and the cloud fractions in the upper and
56    ! lower layers
57    real(jprb), intent(in) :: beta, frac1, frac2
58
59    real(jprb)             :: beta2alpha
60
61    ! Absolute difference in cloud fraction
62    real(jprb)             :: frac_diff
63
64    if (beta < 1.0_jprb) then
65      frac_diff = abs(frac1-frac2)
66      beta2alpha = beta &
67           &  + (1.0_jprb-beta)*frac_diff &
68           &  / (frac_diff + 1.0_jprb/beta - 1.0_jprb)
69    else
70      beta2alpha = 1.0_jprb
71    end if
72
73  end function beta2alpha
74
75
76  !---------------------------------------------------------------------
77  ! Compute total cloud cover according to the specified overlap
78  ! rule. This can be used to compute the high, mid and low cloud
79  ! cover by passing in subsets of the cloud fraction array
80  function cloud_cover(nlev, i_overlap_scheme, frac, overlap_param, &
81       &               is_beta_overlap)
82
83    implicit none
84   
85    ! Number of levels and the overlap scheme to be applied
86    integer, intent(in)    :: nlev, i_overlap_scheme
87
88    ! Cloud fraction and the overlap parameter between adjacent pairs
89    ! of levels
90    real(jprb), intent(in) :: frac(nlev), overlap_param(nlev-1)
91
92    ! Do we use the "beta" overlap scheme of Shonk et al. (2010)?
93    ! Default is false.
94    logical, intent(in), optional :: is_beta_overlap
95
96    ! Return cloud cover
97    real(jprb)             :: cloud_cover
98
99    ! Cumulative cloud cover from TOA to the base of each layer
100    real(jprb) :: cum_cloud_cover(nlev)
101
102    ! Cloud cover of a pair of layers
103    real(jprb) :: pair_cloud_cover(nlev-1)
104
105    if (i_overlap_scheme == IOverlapExponentialRandom) then
106      call cum_cloud_cover_exp_ran(nlev, frac, overlap_param, &
107           &   cum_cloud_cover, pair_cloud_cover, is_beta_overlap)
108    else if (i_overlap_scheme == IOverlapExponential) then
109      call cum_cloud_cover_exp_exp(nlev, frac, overlap_param, &
110           &   cum_cloud_cover, pair_cloud_cover, is_beta_overlap)
111    else
112      call cum_cloud_cover_max_ran(nlev, frac, cum_cloud_cover, &
113           &                       pair_cloud_cover)
114    end if
115
116    cloud_cover = cum_cloud_cover(nlev)
117
118  end function cloud_cover
119
120
121  !---------------------------------------------------------------------
122  ! Maximum-random overlap: Geleyn & Hollingsworth formula
123  subroutine cum_cloud_cover_max_ran(nlev, frac, &
124       & cum_cloud_cover, pair_cloud_cover)
125
126    use yomhook,  only           : lhook, dr_hook, jphook
127
128    implicit none
129
130    ! Inputs
131    integer, intent(in)     :: nlev  ! number of model levels
132
133    ! Cloud fraction on full levels
134    real(jprb), intent(in)  :: frac(nlev)
135
136    ! Outputs
137
138    ! Cumulative cloud cover from TOA to the base of each layer
139    real(jprb), intent(out) :: cum_cloud_cover(nlev)
140
141    ! Cloud cover of a pair of layers
142    real(jprb), intent(out) :: pair_cloud_cover(nlev-1)
143
144    ! Local variables
145
146    ! Cumulative product needed in computation of total_cloud_cover
147    real(jprb) :: cum_product
148
149    ! Loop index for model level
150    integer :: jlev
151
152    real(jphook) :: hook_handle
153
154    if (lhook) call dr_hook('radiation_cloud_cover:cum_cloud_cover_max_ran',0,hook_handle)
155
156    ! Loop to compute total cloud cover and the cumulative cloud cover
157
158    ! down to the base of each layer
159    cum_product = 1.0_jprb - frac(1)
160    cum_cloud_cover(1) = frac(1)
161    do jlev = 1,nlev-1
162      ! Compute the combined cloud cover of layers jlev and jlev+1
163      pair_cloud_cover(jlev) = max(frac(jlev),frac(jlev+1))
164
165      if (frac(jlev) >= MaxCloudFrac) then
166        ! Cloud cover has reached one
167        cum_product = 0.0_jprb
168      else
169        cum_product = cum_product * (1.0_jprb - pair_cloud_cover(jlev)) &
170             &  / (1.0_jprb - frac(jlev))
171      end if
172      cum_cloud_cover(jlev+1) = 1.0_jprb - cum_product
173    end do
174
175    if (lhook) call dr_hook('radiation_cloud_cover:cum_cloud_cover_max_ran',1,hook_handle)
176
177  end subroutine cum_cloud_cover_max_ran
178 
179
180  !---------------------------------------------------------------------
181  ! Exponential-random overlap: exponential overlap for contiguous
182  ! clouds, random overlap for non-contiguous clouds
183  subroutine cum_cloud_cover_exp_ran(nlev, frac, overlap_param, &
184       & cum_cloud_cover, pair_cloud_cover, is_beta_overlap)
185
186    use yomhook,  only           : lhook, dr_hook, jphook
187
188    implicit none
189
190    ! Inputs
191    integer, intent(in)     :: nlev  ! number of model levels
192
193    ! Cloud fraction on full levels
194    real(jprb), intent(in)  :: frac(nlev)
195
196    ! Cloud overlap parameter for interfaces between model layers,
197    ! where 0 indicates random overlap and 1 indicates maximum-random
198    ! overlap
199    real(jprb), intent(in)  :: overlap_param(nlev-1)
200
201    ! This routine has been coded using the "alpha" overlap parameter
202    ! of Hogan and Illingworth (2000). If the following logical is
203    ! present and true then the input is interpretted to be the "beta"
204    ! overlap parameter of Shonk et al. (2010), and needs to be
205    ! converted to alpha.
206    logical, intent(in), optional :: is_beta_overlap
207
208    ! Outputs
209
210    ! Cumulative cloud cover from TOA to the base of each layer
211    real(jprb), intent(out) :: cum_cloud_cover(nlev)
212
213    ! Cloud cover of a pair of layers
214    real(jprb), intent(out) :: pair_cloud_cover(nlev-1)
215
216    ! Local variables
217
218    ! Cumulative product needed in computation of total_cloud_cover
219    real(jprb) :: cum_product
220
221    ! "Alpha" overlap parameter
222    real(jprb) :: overlap_alpha
223    logical    :: do_overlap_conversion
224
225    ! Loop index for model level
226    integer :: jlev
227
228    real(jphook) :: hook_handle
229
230    if (lhook) call dr_hook('radiation_cloud_cover:cum_cloud_cover_exp_ran',0,hook_handle)
231
232   
233    if (present(is_beta_overlap)) then
234      do_overlap_conversion = is_beta_overlap
235    else
236      do_overlap_conversion = .false.
237    end if
238
239    ! Loop to compute total cloud cover and the cumulative cloud cover
240
241    ! down to the base of each layer
242    cum_product = 1.0_jprb - frac(1)
243    cum_cloud_cover(1) = frac(1)
244    do jlev = 1,nlev-1
245      ! Convert to "alpha" overlap parameter if necessary
246      if (do_overlap_conversion) then
247        overlap_alpha = beta2alpha(overlap_param(jlev), &
248             &                     frac(jlev), frac(jlev+1))
249      else
250        overlap_alpha = overlap_param(jlev)
251      end if
252
253      ! Compute the combined cloud cover of layers jlev and jlev+1
254      pair_cloud_cover(jlev) = overlap_alpha*max(frac(jlev),frac(jlev+1)) &
255           &  + (1.0_jprb - overlap_alpha) &
256           &  * (frac(jlev)+frac(jlev+1)-frac(jlev)*frac(jlev+1))
257! Added for DWD (2020)
258#ifdef DWD_VECTOR_OPTIMIZATIONS
259    end do
260    do jlev = 1,nlev-1
261#endif
262      if (frac(jlev) >= MaxCloudFrac) then
263        ! Cloud cover has reached one
264        cum_product = 0.0_jprb
265      else
266        cum_product = cum_product * (1.0_jprb - pair_cloud_cover(jlev)) &
267             &  / (1.0_jprb - frac(jlev))
268      end if
269      cum_cloud_cover(jlev+1) = 1.0_jprb - cum_product
270    end do
271
272
273    if (lhook) call dr_hook('radiation_cloud_cover:cum_cloud_cover_exp_ran',1,hook_handle)
274
275  end subroutine cum_cloud_cover_exp_ran
276
277
278
279  !---------------------------------------------------------------------
280  ! Exponential-exponential overlap: exponential overlap for both
281  ! contiguous and non-contiguous clouds. This is the result of the
282  ! simple Raisanen cloud generator, but unfortunately it has no
283  ! (known) analytic formula for the total cloud cover, or the
284  ! cumulative cloud cover.  In partially cloudy columns, The McICA
285  ! scheme needs this info in order to devote all the cloudy g-points
286  ! to columns containing cloud, which reduces McICA noise. The
287  ! following routine provides an approximate estimate of cumulative
288  ! cloud cover consistent with the exponential-exponential scheme.
289  subroutine cum_cloud_cover_exp_exp(nlev, frac, overlap_param, &
290       & cum_cloud_cover, pair_cloud_cover, is_beta_overlap)
291
292    use yomhook,  only           : lhook, dr_hook, jphook
293
294    implicit none
295
296    ! Inputs
297    integer, intent(in)     :: nlev  ! number of model levels
298
299    ! Cloud fraction on full levels
300    real(jprb), intent(in)  :: frac(nlev)
301
302    ! Cloud overlap parameter for interfaces between model layers,
303    ! where 0 indicates random overlap and 1 indicates maximum-random
304    ! overlap
305    real(jprb), intent(in)  :: overlap_param(nlev-1)
306
307    ! This routine has been coded using the "alpha" overlap parameter
308    ! of Hogan and Illingworth (2000). If the following logical is
309    ! present and true then the input is interpretted to be the "beta"
310    ! overlap parameter of Shonk et al. (2010), and needs to be
311    ! converted to alpha.
312    logical, intent(in), optional :: is_beta_overlap
313
314    ! Outputs
315
316    ! Cumulative cloud cover from TOA to the base of each layer
317    real(jprb), intent(out) :: cum_cloud_cover(nlev)
318
319    ! Cloud cover of a pair of layers
320    real(jprb), intent(out) :: pair_cloud_cover(nlev-1)
321
322    ! Local variables
323
324    ! If this routine is called from the radiation_interface tree then
325    ! very low cloud fractions have already been set to zero, but if
326    ! it is called as a cloud cover diagnostic then this can't be
327    ! guaranteed so a small non-zero numbers is required
328    real(jprb), parameter :: min_frac = 1.0e-6_jprb
329
330    ! "Alpha" overlap parameter
331    real(jprb) :: overlap_alpha(nlev-1)
332    logical    :: do_overlap_conversion
333
334    ! Variables describing "concave cloud objects", i.e. contiguous
335    ! layers where the cloud fraction monotonically increases then
336    ! monotonically decreases
337
338    ! Number of objects
339    integer    :: nobj
340
341    ! Indices to the location of the top, maximum cloud fraction, and
342    ! base, of each concave cloud object
343    integer    :: i_top_obj(nlev)
344    integer    :: i_max_obj(nlev)
345    integer    :: i_base_obj(nlev)
346    ! Poor-man's linked list to allow for deletion of objects: this
347    ! variable points to the index of the next active object
348    integer    :: i_next_obj(nlev)
349
350    ! Cloud cover of object
351    real(jprb) :: cc_obj(nlev)
352
353    ! Overlap parameter between objects
354    real(jprb) :: alpha_obj(nlev)
355
356    ! Do (while) loop index for model level
357    integer :: jlev
358
359    ! Do loop index for object
360    integer :: jobj
361
362    ! Maximum correlation between adjacent objects
363    real(jprb) :: alpha_max
364
365    ! Indices to pair of objects
366    integer    :: iobj1, iobj2
367
368    ! Combined cloud cover of pair of objects, and scaling to modify
369    ! cumulative cloud cover of lower layer
370    real(jprb) :: cc_pair, scaling
371
372    real(jphook) :: hook_handle
373
374    if (lhook) call dr_hook('radiation_cloud_cover:cum_cloud_cover_exp_exp',0,hook_handle)
375
376   
377    if (present(is_beta_overlap)) then
378      do_overlap_conversion = is_beta_overlap
379    else
380      do_overlap_conversion = .false.
381    end if
382
383    ! Loop down through atmosphere to locate objects and compute their
384    ! basic properties
385    jlev = 1
386    nobj = 0
387    do while (jlev <= nlev)
388      if (frac(jlev) > min_frac) then
389        ! Starting a new object: store its top
390        nobj = nobj + 1
391        i_top_obj(nobj) = jlev;
392        ! Find its maximum cloud fraction
393        jlev = jlev + 1
394        do while (jlev <= nlev)
395          if (frac(jlev) < frac(jlev-1)) then
396            exit
397          end if
398          jlev = jlev + 1
399        end do
400        i_max_obj(nobj) = jlev - 1
401        ! Find its base
402        do while (jlev <= nlev)
403          if (frac(jlev) > frac(jlev-1) .or. frac(jlev) <= min_frac) then
404            exit
405          end if
406          jlev = jlev + 1
407        end do
408        ! In the case of cloud fraction profile starting from the top
409        ! like this: 0.1 0.2 0.1 0.2 0.1, we may want the object grouping
410        ! to be (0.1 0.2) (0.1 0.2 0.1), not (0.1 0.2 0.1) (0.2 0.1),
411        ! in which case the following should be uncommented
412        !        if (jlev < nlev) then
413        !          if (frac(jlev) > frac(jlev-1)) then
414        !            jlev = jlev - 1
415        !          end if
416        !        end if
417        i_base_obj(nobj) = jlev - 1
418        ! Index to the next object
419        i_next_obj(nobj) = nobj+1
420      else
421        jlev = jlev + 1
422      end if
423    end do
424
425    ! Array assignments
426    cum_cloud_cover = 0.0_jprb
427    pair_cloud_cover = 0.0_jprb
428
429    if (nobj > 0) then
430      ! Only do any more work if there is cloud present
431     
432      ! To minimize the potential calls to beta2alpha, we do all the
433      ! computations related to overlap parameter here
434      if (.not. do_overlap_conversion) then
435
436        ! Compute the combined cloud cover of pairs of layers
437        do jlev = 1,nlev-1
438          pair_cloud_cover(jlev) &
439               &  = overlap_param(jlev)*max(frac(jlev),frac(jlev+1)) &
440               &  + (1.0_jprb - overlap_param(jlev)) &
441               &  * (frac(jlev)+frac(jlev+1)-frac(jlev)*frac(jlev+1))
442        end do
443        ! Estimate the effective overlap parameter "alpha_obj" between
444        ! adjacent objects as the product of the layerwise overlap
445        ! parameters between their layers of maximum cloud fraction
446        do jobj = 1,nobj-1
447          alpha_obj(jobj) &
448               &  = product(overlap_param(i_max_obj(jobj):i_max_obj(jobj+1)-1))
449        end do
450
451      else
452
453        ! Convert Shonk et al overlap parameter to Hogan and
454        ! Illingworth definition
455        overlap_alpha = beta2alpha(overlap_param, &
456             &                     frac(1:nlev-1), frac(2:nlev))
457        ! Compute the combined cloud cover of pairs of layers
458        do jlev = 1,nlev-1
459          pair_cloud_cover(jlev) &
460               &  = overlap_alpha(jlev)*max(frac(jlev),frac(jlev+1)) &
461               &  + (1.0_jprb - overlap_alpha(jlev)) &
462               &  * (frac(jlev)+frac(jlev+1)-frac(jlev)*frac(jlev+1))         
463        end do
464        ! Estimate the effective overlap parameter "alpha_obj" between
465        ! adjacent objects as the product of the layerwise overlap
466        ! parameters between their layers of maximum cloud fraction
467        do jobj = 1,nobj-1
468          alpha_obj(jobj) &
469               &  = product(overlap_alpha(i_max_obj(jobj):i_max_obj(jobj+1)-1))
470        end do
471
472      end if
473
474
475      ! Compute the cumulative cloud cover working down from the top
476      ! of each object: this will later be converted to the cumulative
477      ! cloud cover working down from TOA
478      do jobj = 1,nobj
479        cum_cloud_cover(i_top_obj(jobj)) = frac(i_top_obj(jobj))
480        do jlev = i_top_obj(jobj), i_base_obj(jobj)-1
481          if (frac(jlev) >= MaxCloudFrac) then
482            ! Cloud cover has reached one
483            cum_cloud_cover(jlev+1) = 1.0_jprb
484          else
485            cum_cloud_cover(jlev+1) = 1.0_jprb &
486                 &  - (1.0_jprb - cum_cloud_cover(jlev)) &
487                 &  * (1.0_jprb - pair_cloud_cover(jlev))  &
488                 &  / (1.0_jprb - frac(jlev))
489          end if
490        end do
491        cc_obj(jobj) = cum_cloud_cover(i_base_obj(jobj))
492      end do
493
494      iobj1 = 1
495
496      ! Sequentially combine objects until there is only one left
497      ! covering the full vertical extent of clouds in the profile
498      do while (nobj > 1)
499        ! Find the most correlated adjacent pair of objects
500        alpha_max = 0.0_jprb
501
502        ! Need to re-initialize iobj1 here in case alpha_obj(:)==0.0,
503        ! which would mean that the "if" statement in the following
504        ! loop would never get triggered
505        iobj1 = 1
506
507        jobj = 1
508        do while (jobj < nobj)
509          if (alpha_obj(jobj) > alpha_max) then
510            alpha_max = alpha_obj(jobj)
511            iobj1 = jobj
512          end if
513          jobj = i_next_obj(jobj)
514        end do
515
516        ! iobj1 is the index to the first object in the pair, set
517        ! iobj2 to the second
518        iobj2 = i_next_obj(iobj1)
519
520        ! Set the cumulative cloud cover in the clear-sky gap between
521        ! the objects to the value at the base of the upper object
522        cum_cloud_cover(i_base_obj(iobj1)+1:i_top_obj(iobj2)-1) &
523             &  = cum_cloud_cover(i_base_obj(iobj1))
524
525        ! Calculate the combined cloud cover of the pair of objects
526        cc_pair = alpha_obj(iobj1)*max(cc_obj(iobj1), cc_obj(iobj2)) &
527             &  + (1.0_jprb - alpha_obj(iobj1)) &
528             &  * (cc_obj(iobj1) + cc_obj(iobj2) - cc_obj(iobj1)*cc_obj(iobj2))
529        scaling = min(max((cc_pair-cc_obj(iobj1)) / max(min_frac, cc_obj(iobj2)), &
530             &            0.0_jprb), &
531             &        1.0_jprb)
532       
533        ! Scale the combined cloud cover of the lower object to
534        ! account for its overlap with the upper object
535        do jlev = i_top_obj(iobj2),i_base_obj(iobj2)
536          cum_cloud_cover(jlev) = cum_cloud_cover(i_base_obj(iobj1)) &
537               +  cum_cloud_cover(jlev) * scaling
538        end do
539       
540        ! Merge the objects by setting the properties of the upper
541        ! object to the combined properties of both.  Note that
542        ! i_max_obj is not modified because it is no longer needed.
543        cc_obj(iobj1) = cc_pair
544        i_base_obj(iobj1) = i_base_obj(iobj2)
545        i_next_obj(iobj1) = i_next_obj(iobj2)
546        alpha_obj(iobj1)  = alpha_obj(iobj2)
547        nobj = nobj - 1
548      end do
549
550      ! Finish off the total cloud cover below cloud
551      cum_cloud_cover(i_base_obj(iobj1)+1:nlev) &
552           &  = cum_cloud_cover(i_base_obj(iobj1))
553
554      ! Ensure that the combined cloud cover of pairs of layers is
555      ! consistent with the overhang
556      do jlev = 1,nlev-1
557        pair_cloud_cover(jlev) = max(pair_cloud_cover(jlev), &
558             &     frac(jlev)+cum_cloud_cover(jlev+1)-cum_cloud_cover(jlev))
559      end do
560
561      ! Sometimes round-off error can lead to cloud cover just above
562      ! one, which in turn can lead to direct shortwave fluxes just
563      ! below zero
564      cum_cloud_cover = min(cum_cloud_cover, 1.0_jprb)
565
566    end if ! cloud is present in profile
567
568    if (lhook) call dr_hook('radiation_cloud_cover:cum_cloud_cover_exp_exp',1,hook_handle)
569
570  end subroutine cum_cloud_cover_exp_exp
571
572end module radiation_cloud_cover
Note: See TracBrowser for help on using the repository browser.