source: LMDZ6/trunk/libf/phylmd/ecrad/radiation_cloud_generator.F90 @ 4163

Last change on this file since 4163 was 3908, checked in by idelkadi, 4 years ago

Online implementation of the radiative transfer code ECRAD in the LMDZ model.

  • Inclusion of the ecrad directory containing the sources of the ECRAD code
    • interface routine : radiation_scheme.F90
  • Adaptation of compilation scripts :
    • compilation under CPP key CPP_ECRAD
    • compilation with option "-rad ecard" or "-ecard true"
    • The "-rad old/rtm/ecran" build option will need to replace the "-rrtm true" and "-ecrad true" options in the future.
  • Runing LMDZ simulations with ecrad, you need :
    • logical key iflag_rrtm = 2 in physiq.def
    • namelist_ecrad (DefLists?)
    • the directory "data" containing the configuration files is temporarily placed in ../libfphylmd/ecrad/
  • Compilation and execution are tested in the 1D case. The repository under svn would allow to continue the implementation work: tests, verification of the results, ...
File size: 17.3 KB
Line 
1! radiation_cloud_generator.F90 - Generate water-content or optical-depth scalings for McICA
2!
3! (C) Copyright 2015- 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 clouds for McICA using a method modified from Raisanen et
16! al. (2002)
17!
18! Modifications
19!   2018-02-22  R. Hogan  Call masked version of PDF sampler for speed
20
21module radiation_cloud_generator
22
23  public
24
25contains
26
27  !---------------------------------------------------------------------
28  ! Generate scaling factors for the cloud optical depth to represent
29  ! cloud overlap, the overlap of internal cloud inhomogeneities and
30  ! the fractional standard deviation of these inhomogeneities, for
31  ! use in a Monte Carlo Independent Column Approximation radiation
32  ! scheme. All returned profiles contain cloud, and the total cloud
33  ! cover is also returned, so the calling function can then do a
34  ! weighted average of clear and cloudy skies; this is a way to
35  ! reduce the Monte Carlo noise in profiles with low cloud cover.
36  subroutine cloud_generator(ng, nlev, i_overlap_scheme, &
37       &  iseed, frac_threshold, &
38       &  frac, overlap_param, decorrelation_scaling, &
39       &  fractional_std, pdf_sampler, &
40       &  od_scaling, total_cloud_cover, &
41       &  is_beta_overlap)
42
43    use parkind1, only           : jprb
44    use yomhook,  only           : lhook, dr_hook
45    use radiation_io,   only     : nulerr, radiation_abort
46    use random_numbers_mix, only : randomnumberstream, &
47         initialize_random_numbers, uniform_distribution
48    use radiation_pdf_sampler, only : pdf_sampler_type
49    use radiation_cloud_cover, only : cum_cloud_cover_exp_ran, &
50         &       cum_cloud_cover_max_ran, cum_cloud_cover_exp_exp, &
51         &       IOverlapMaximumRandom, IOverlapExponentialRandom, &
52         &       IOverlapExponential
53
54    implicit none
55
56    ! Inputs
57    integer, intent(in)     :: ng    ! number of g points
58    integer, intent(in)     :: nlev  ! number of model levels
59    integer, intent(in)     :: i_overlap_scheme
60    integer, intent(in)     :: iseed ! seed for random number generator
61
62    ! Only cloud fractions above this threshold are considered to be
63    ! clouds
64    real(jprb), intent(in)  :: frac_threshold
65
66    ! Cloud fraction on full levels
67    real(jprb), intent(in)  :: frac(nlev)
68
69    ! Cloud overlap parameter for interfaces between model layers,
70    ! where 0 indicates random overlap and 1 indicates maximum-random
71    ! overlap
72    real(jprb), intent(in)  :: overlap_param(nlev-1)
73
74    ! Overlap parameter for internal inhomogeneities
75    real(jprb), intent(in)  :: decorrelation_scaling
76
77    ! Fractional standard deviation at each layer
78    real(jprb), intent(in)  :: fractional_std(nlev)
79
80    ! Object for sampling from a lognormal or gamma distribution
81    type(pdf_sampler_type), intent(in) :: pdf_sampler
82
83    ! This routine has been coded using the "alpha" overlap parameter
84    ! of Hogan and Illingworth (2000). If the following logical is
85    ! present and true then the input is interpretted to be the "beta"
86    ! overlap parameter of Shonk et al. (2010), and needs to be
87    ! converted to alpha.
88    logical, intent(in), optional :: is_beta_overlap
89
90    ! Outputs
91
92    ! Cloud optical depth scaling factor, with 0 indicating clear sky
93    real(jprb), intent(out) :: od_scaling(ng,nlev)
94
95    ! Total cloud cover using cloud fraction and overlap parameter
96    real(jprb), intent(out) :: total_cloud_cover
97
98    ! Local variables
99
100    ! Cumulative cloud cover from TOA to the base of each layer
101    real(jprb) :: cum_cloud_cover(nlev)
102
103    ! Scaled random number for finding cloud
104    real(jprb) :: trigger
105
106    ! Uniform deviates between 0 and 1
107    real(jprb) :: rand_top(ng)
108
109    ! Overlap parameter of inhomogeneities
110    real(jprb) :: overlap_param_inhom(nlev-1)
111
112    ! Seed for random number generator and stream for producing random
113    ! numbers
114    type(randomnumberstream) :: random_stream
115   
116    ! First and last cloudy layers
117    integer :: ibegin, iend
118
119    integer :: itrigger
120
121    ! Loop index for model level and g-point
122    integer :: jlev, jg
123
124    ! Cloud cover of a pair of layers, and amount by which cloud at
125    ! next level increases total cloud cover as seen from above
126    real(jprb), dimension(nlev-1) :: pair_cloud_cover, overhang
127
128    real(jprb) :: hook_handle
129
130    if (lhook) call dr_hook('radiation_cloud_generator:cloud_generator',0,hook_handle)
131
132    if (i_overlap_scheme == IOverlapExponentialRandom) then
133      call cum_cloud_cover_exp_ran(nlev, frac, overlap_param, &
134           &   cum_cloud_cover, pair_cloud_cover, is_beta_overlap)
135    else if (i_overlap_scheme == IOverlapMaximumRandom) then
136      call cum_cloud_cover_max_ran(nlev, frac, &
137           &   cum_cloud_cover, pair_cloud_cover)
138    else if (i_overlap_scheme == IOverlapExponential) then
139      call cum_cloud_cover_exp_exp(nlev, frac, overlap_param, &
140           &   cum_cloud_cover, pair_cloud_cover, is_beta_overlap)
141    else
142      write(nulerr,'(a)') '*** Error: cloud overlap scheme not recognised'
143      call radiation_abort()
144    end if
145
146    total_cloud_cover = cum_cloud_cover(nlev);
147    do jlev = 1,nlev-1
148      overhang(jlev) = cum_cloud_cover(jlev+1)-cum_cloud_cover(jlev)
149    end do
150
151    if (total_cloud_cover < frac_threshold) then
152      ! Treat column as clear sky: calling function therefore will not
153      ! use od_scaling so we don't need to calculate it
154      total_cloud_cover = 0.0_jprb
155
156    else
157      ! Cloud is present: need to calculate od_scaling
158
159      ! Find range of cloudy layers
160      jlev = 1
161      do while (frac(jlev) <= 0.0_jprb)
162        jlev = jlev + 1
163      end do
164      ibegin = jlev
165      iend = jlev
166      do jlev = jlev+1,nlev
167        if (frac(jlev) > 0.0_jprb) then
168          iend = jlev
169        end if
170      end do
171
172      ! Set overlap parameter of inhomogeneities
173      overlap_param_inhom = overlap_param
174
175      do jlev = ibegin,iend-1
176        if (overlap_param(jlev) > 0.0_jprb) then
177          overlap_param_inhom(jlev) &
178               &  = overlap_param(jlev)**(1.0_jprb/decorrelation_scaling)
179        end if
180      end do
181
182      ! Reset optical depth scaling to clear skies
183      od_scaling = 0.0_jprb
184
185      ! Expensive operation: initialize random number generator for
186      ! this column
187      call initialize_random_numbers(iseed, random_stream)
188
189      ! Compute ng random numbers to use to locate cloud top
190      call uniform_distribution(rand_top, random_stream)
191
192      ! Loop over ng columns
193      do jg = 1,ng
194        ! Find the cloud top height corresponding to the current
195        ! random number, and store in itrigger
196        trigger = rand_top(jg) * total_cloud_cover
197        jlev = ibegin
198        do while (trigger > cum_cloud_cover(jlev) .and. jlev < iend)
199          jlev = jlev + 1
200        end do
201        itrigger = jlev
202
203        if (i_overlap_scheme /= IOverlapExponential) then
204          call generate_column_exp_ran(ng, nlev, jg, random_stream, pdf_sampler, &
205               &  frac, pair_cloud_cover, &
206               &  cum_cloud_cover, overhang, fractional_std, overlap_param_inhom, &
207               &  itrigger, iend, od_scaling)
208        else
209          call generate_column_exp_exp(ng, nlev, jg, random_stream, pdf_sampler, &
210               &  frac, pair_cloud_cover, &
211               &  cum_cloud_cover, overhang, fractional_std, overlap_param_inhom, &
212               &  itrigger, iend, od_scaling)
213        end if
214       
215      end do
216
217    end if
218
219
220    if (lhook) call dr_hook('radiation_cloud_generator:cloud_generator',1,hook_handle)
221
222  end subroutine cloud_generator
223
224
225  !---------------------------------------------------------------------
226  ! Generate a column of optical depth scalings using
227  ! exponential-random overlap (which includes maximum-random overlap
228  ! as a limiting case)
229  subroutine generate_column_exp_ran(ng, nlev, ig, random_stream, pdf_sampler, &
230       &  frac, pair_cloud_cover, &
231       &  cum_cloud_cover, overhang, fractional_std, overlap_param_inhom, &
232       &  itrigger, iend, od_scaling)
233
234    use parkind1,              only : jprb
235    use radiation_pdf_sampler, only : pdf_sampler_type
236    use random_numbers_mix,    only : randomnumberstream, &
237         initialize_random_numbers, uniform_distribution
238
239
240    implicit none
241
242    ! Number of g points / columns, and number of current column
243    integer, intent(in) :: ng, ig
244
245    ! Number of levels
246    integer, intent(in) :: nlev
247
248    ! Stream for producing random numbers
249    type(randomnumberstream), intent(inout) :: random_stream
250
251    ! Object for sampling from a lognormal or gamma distribution
252    type(pdf_sampler_type), intent(in) :: pdf_sampler
253
254    ! Cloud fraction, cumulative cloud cover and fractional standard
255    ! deviation in each layer
256    real(jprb), intent(in), dimension(nlev) :: frac, cum_cloud_cover, fractional_std
257
258    ! Cloud cover of a pair of layers, and amount by which cloud at
259    ! next level increases total cloud cover as seen from above
260    real(jprb), intent(in), dimension(nlev-1) :: pair_cloud_cover, overhang
261
262    ! Overlap parameter of inhomogeneities
263    real(jprb), intent(in), dimension(nlev-1) :: overlap_param_inhom
264
265    ! Top of highest cloudy layer (in this subcolumn) and base of
266    ! lowest
267    integer, intent(in) :: itrigger, iend
268
269    ! Optical depth scaling to output
270    real(jprb), intent(inout), dimension(ng,nlev) :: od_scaling
271
272    ! Height indices
273    integer :: jlev, jcloud
274
275    ! Number of contiguous cloudy layers for which to compute optical
276    ! depth scaling
277    integer :: n_layers_to_scale
278
279    integer :: iy
280
281    ! Is it time to fill the od_scaling variable?
282    logical :: do_fill_od_scaling
283
284    real(jprb) :: rand_cloud(nlev)
285    real(jprb) :: rand_inhom1(nlev), rand_inhom2(nlev)
286
287    ! So far our vertically contiguous cloud contains only one layer
288    n_layers_to_scale = 1
289    iy = 0
290
291    ! Locate the clouds below this layer: first generate some more
292    ! random numbers
293    call uniform_distribution(rand_cloud(1:(iend+1-itrigger)),random_stream)
294
295    ! Loop from the layer below the local cloud top down to the
296    ! bottom-most cloudy layer
297    do jlev = itrigger+1,iend+1
298      do_fill_od_scaling = .false.
299      if (jlev <= iend) then
300        iy = iy+1
301        if (n_layers_to_scale > 0) then
302          ! There is a cloud above, in which case the probability
303          ! of cloud in the layer below is as follows
304          if (rand_cloud(iy)*frac(jlev-1) &
305               &  < frac(jlev) + frac(jlev-1) - pair_cloud_cover(jlev-1)) then
306            ! Add another cloudy layer
307            n_layers_to_scale = n_layers_to_scale + 1
308          else
309            ! Reached the end of a contiguous set of cloudy layers and
310            ! will compute the optical depth scaling immediately.
311            do_fill_od_scaling = .true.
312          end if
313        else
314          ! There is clear-sky above, in which case the
315          ! probability of cloud in the layer below is as follows
316          if (rand_cloud(iy)*(cum_cloud_cover(jlev-1) - frac(jlev-1)) &
317               &  < pair_cloud_cover(jlev-1) - overhang(jlev-1) - frac(jlev-1)) then
318            ! A new cloud top
319            n_layers_to_scale = 1
320          end if
321        end if
322      else
323        ! We are at the bottom of the cloudy layers in the model,
324        ! so in a moment need to populate the od_scaling array
325        do_fill_od_scaling = .true.
326      end if
327
328      if (do_fill_od_scaling) then
329        ! We have a contiguous range of layers for which we
330        ! compute the od_scaling elements using some random
331        ! numbers
332        call uniform_distribution(rand_inhom1(1:n_layers_to_scale),random_stream)
333        call uniform_distribution(rand_inhom2(1:n_layers_to_scale),random_stream)
334
335        ! Loop through the sequence of cloudy layers
336        do jcloud = 2,n_layers_to_scale
337          ! Use second random number, and inhomogeneity overlap
338          ! parameter, to decide whether the first random number
339          ! should be repeated (corresponding to maximum overlap)
340          ! or not (corresponding to random overlap)
341          if (rand_inhom2(jcloud) &
342               &  < overlap_param_inhom(jlev-n_layers_to_scale+jcloud-2)) then
343            rand_inhom1(jcloud) = rand_inhom1(jcloud-1)
344          end if
345        end do
346       
347        ! Sample from a lognormal or gamma distribution to obtain
348        ! the optical depth scalings
349        call pdf_sampler%sample(fractional_std(jlev-n_layers_to_scale:jlev-1), &
350             & rand_inhom1(1:n_layers_to_scale), od_scaling(ig,jlev-n_layers_to_scale:jlev-1))
351
352        n_layers_to_scale = 0
353      end if
354         
355    end do
356
357  end subroutine generate_column_exp_ran
358
359
360
361  !---------------------------------------------------------------------
362  ! Generate a column of optical depth scalings using
363  ! exponential-exponential overlap
364  subroutine generate_column_exp_exp(ng, nlev, ig, random_stream, pdf_sampler, &
365       &  frac, pair_cloud_cover, &
366       &  cum_cloud_cover, overhang, fractional_std, overlap_param_inhom, &
367       &  itrigger, iend, od_scaling)
368
369    use parkind1,              only : jprb
370    use radiation_pdf_sampler, only : pdf_sampler_type
371    use random_numbers_mix,    only : randomnumberstream, &
372         initialize_random_numbers, uniform_distribution
373
374    implicit none
375
376    ! Number of g points / columns, and number of current column
377    integer, intent(in) :: ng, ig
378
379    ! Number of levels
380    integer, intent(in) :: nlev
381
382    ! Stream for producing random numbers
383    type(randomnumberstream), intent(inout) :: random_stream
384
385    ! Object for sampling from a lognormal or gamma distribution
386    type(pdf_sampler_type), intent(in) :: pdf_sampler
387
388    ! Cloud fraction, cumulative cloud cover and fractional standard
389    ! deviation in each layer
390    real(jprb), intent(in), dimension(nlev) :: frac, cum_cloud_cover, fractional_std
391
392    ! Cloud cover of a pair of layers, and amount by which cloud at
393    ! next level increases total cloud cover as seen from above
394    real(jprb), intent(in), dimension(nlev-1) :: pair_cloud_cover, overhang
395
396    ! Overlap parameter of inhomogeneities
397    real(jprb), intent(in), dimension(nlev-1) :: overlap_param_inhom
398
399    ! Top of highest cloudy layer (in this subcolumn) and base of
400    ! lowest
401    integer, intent(in) :: itrigger, iend
402
403    ! Optical depth scaling to output
404    real(jprb), intent(inout), dimension(ng,nlev) :: od_scaling
405
406    ! Height indices
407    integer :: jlev, jcloud
408
409    integer :: iy
410
411    real(jprb) :: rand_cloud(nlev)
412    real(jprb) :: rand_inhom1(nlev), rand_inhom2(nlev)
413
414    ! For each column analysed, this vector locates the clouds. It is
415    ! only actually used for Exp-Exp overlap
416    logical :: is_cloudy(nlev)
417
418    ! Number of contiguous cloudy layers for which to compute optical
419    ! depth scaling
420    integer :: n_layers_to_scale
421
422    iy = 0
423
424    is_cloudy = .false.
425    is_cloudy(itrigger) = .true.
426
427    ! Locate the clouds below this layer: first generate some more
428    ! random numbers
429    call uniform_distribution(rand_cloud(1:(iend+1-itrigger)),random_stream)
430
431    ! Loop from the layer below the local cloud top down to the
432    ! bottom-most cloudy layer
433    do jlev = itrigger+1,iend
434      iy = iy+1
435      if (is_cloudy(jlev-1)) then
436        ! There is a cloud above, in which case the probability
437        ! of cloud in the layer below is as follows
438        if (rand_cloud(iy)*frac(jlev-1) &
439             &  < frac(jlev) + frac(jlev-1) - pair_cloud_cover(jlev-1)) then
440          ! Add another cloudy layer
441          is_cloudy(jlev) = .true.
442        end if
443      else
444        ! There is clear-sky above, in which case the
445        ! probability of cloud in the layer below is as follows
446        if (rand_cloud(iy)*(cum_cloud_cover(jlev-1) - frac(jlev-1)) &
447             &  < pair_cloud_cover(jlev-1) - overhang(jlev-1) - frac(jlev-1)) then
448            ! A new cloud top
449          is_cloudy(jlev) = .true.
450        end if
451      end if
452    end do
453
454    ! We have a contiguous range of layers for which we compute the
455    ! od_scaling elements using some random numbers
456
457    ! In the Exp-Exp overlap scheme we do all layers at once
458    n_layers_to_scale = iend+1 - itrigger
459       
460    call uniform_distribution(rand_inhom1(1:n_layers_to_scale),random_stream)
461    call uniform_distribution(rand_inhom2(1:n_layers_to_scale),random_stream)
462       
463    ! Loop through the sequence of cloudy layers
464    do jcloud = 2,n_layers_to_scale
465      ! Use second random number, and inhomogeneity overlap
466      ! parameter, to decide whether the first random number
467      ! should be repeated (corresponding to maximum overlap)
468      ! or not (corresponding to random overlap)
469      if (rand_inhom2(jcloud) &
470           &  < overlap_param_inhom(iend-n_layers_to_scale+jcloud-1)) then
471        rand_inhom1(jcloud) = rand_inhom1(jcloud-1)
472      end if
473    end do
474       
475    ! Sample from a lognormal or gamma distribution to obtain the
476    ! optical depth scalings, calling the faster masked version and
477    ! assuming values outside the range itrigger:iend are already zero
478    call pdf_sampler%masked_sample(n_layers_to_scale, &
479         &  fractional_std(itrigger:iend), &
480         &  rand_inhom1(1:n_layers_to_scale), od_scaling(ig,itrigger:iend), &
481         &  is_cloudy(itrigger:iend))
482       
483  end subroutine generate_column_exp_exp
484
485end module radiation_cloud_generator
Note: See TracBrowser for help on using the repository browser.