source: LMDZ6/trunk/libf/phylmd/ecrad.v1.5.1/radiation_cloud_cover.F90 @ 4894

Last change on this file since 4894 was 4489, checked in by lguez, 20 months ago

Merge LMDZ_ECRad branch back into trunk!

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