source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_two_stream.F90 @ 4773

Last change on this file since 4773 was 4773, checked in by idelkadi, 5 months ago
  • Update of Ecrad in LMDZ The same organization of the Ecrad offline version is retained in order to facilitate the updating of Ecrad in LMDZ and the comparison between online and offline results. version 1.6.1 of Ecrad (https://github.com/lguez/ecrad.git)
  • Implementation of the double call of Ecrad in LMDZ


File size: 31.2 KB
Line 
1! radiation_two_stream.F90 - Compute two-stream coefficients
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!   2017-05-04  P Dueben/R Hogan  Use JPRD where double precision essential
17!   2017-07-12  R Hogan  Optimized LW coeffs in low optical depth case
18!   2017-07-26  R Hogan  Added calc_frac_scattered_diffuse_sw routine
19!   2017-10-23  R Hogan  Renamed single-character variables
20!   2021-02-19  R Hogan  Security for shortwave singularity
21!   2022-11-22  P Ukkonen/R Hogan  Single precision uses no double precision
22!   2023-09-28  R Hogan  Increased security for single-precision SW "k"
23
24module radiation_two_stream
25
26  use parkind1, only : jprb, jprd
27
28  implicit none
29  public
30
31  ! Elsasser's factor: the effective factor by which the zenith
32  ! optical depth needs to be multiplied to account for longwave
33  ! transmission at all angles through the atmosphere.  Alternatively
34  ! think of acos(1/lw_diffusivity) to be the effective zenith angle
35  ! of longwave radiation.
36  real(jprd), parameter :: LwDiffusivity   = 1.66_jprd
37  real(jprb), parameter :: LwDiffusivityWP = 1.66_jprb ! Working precision version
38
39  ! The routines in this module can be called millions of times, so
40  ! calling Dr Hook for each one may be a significant overhead.
41  ! Uncomment the following to turn Dr Hook on.
42!#define DO_DR_HOOK_TWO_STREAM
43
44contains
45
46#ifdef FAST_EXPONENTIAL
47  !---------------------------------------------------------------------
48  ! Fast exponential for negative arguments: a Pade approximant that
49  ! doesn't go negative for negative arguments, applied to arg/8, and
50  ! the result is then squared three times
51  elemental function exp_fast(arg) result(ex)
52    real(jprd) :: arg, ex
53    ex = 1.0_jprd / (1.0_jprd + arg*(-0.125_jprd &
54         + arg*(0.0078125_jprd - 0.000325520833333333_jprd * arg)))
55    ex = ex*ex
56    ex = ex*ex
57    ex = ex*ex
58  end function exp_fast
59#else
60#define exp_fast exp
61#endif
62
63  !---------------------------------------------------------------------
64  ! Calculate the two-stream coefficients gamma1 and gamma2 for the
65  ! longwave
66  subroutine calc_two_stream_gammas_lw(ng, ssa, g, &
67       &                               gamma1, gamma2)
68
69#ifdef DO_DR_HOOK_TWO_STREAM
70    use yomhook, only : lhook, dr_hook, jphook
71#endif
72
73    integer, intent(in) :: ng
74    ! Sngle scattering albedo and asymmetry factor:
75    real(jprb), intent(in),  dimension(:) :: ssa, g
76    real(jprb), intent(out), dimension(:) :: gamma1, gamma2
77
78    real(jprb) :: factor
79
80    integer    :: jg
81
82#ifdef DO_DR_HOOK_TWO_STREAM
83    real(jphook) :: hook_handle
84
85    if (lhook) call dr_hook('radiation_two_stream:calc_two_stream_gammas_lw',0,hook_handle)
86#endif
87! Added for DWD (2020)
88!NEC$ shortloop
89    do jg = 1, ng
90      ! Fu et al. (1997), Eq 2.9 and 2.10:
91      !      gamma1(jg) = LwDiffusivity * (1.0_jprb - 0.5_jprb*ssa(jg) &
92      !           &                    * (1.0_jprb + g(jg)))
93      !      gamma2(jg) = LwDiffusivity * 0.5_jprb * ssa(jg) &
94      !           &                    * (1.0_jprb - g(jg))
95      ! Reduce number of multiplications
96      factor = (LwDiffusivity * 0.5_jprb) * ssa(jg)
97      gamma1(jg) = LwDiffusivity - factor*(1.0_jprb + g(jg))
98      gamma2(jg) = factor * (1.0_jprb - g(jg))
99    end do
100
101#ifdef DO_DR_HOOK_TWO_STREAM
102    if (lhook) call dr_hook('radiation_two_stream:calc_two_stream_gammas_lw',1,hook_handle)
103#endif
104
105  end subroutine calc_two_stream_gammas_lw
106
107
108  !---------------------------------------------------------------------
109  ! Calculate the two-stream coefficients gamma1-gamma4 in the
110  ! shortwave
111  subroutine calc_two_stream_gammas_sw(ng, mu0, ssa, g, &
112       &                               gamma1, gamma2, gamma3)
113
114#ifdef DO_DR_HOOK_TWO_STREAM
115    use yomhook, only : lhook, dr_hook, jphook
116#endif
117
118    integer, intent(in) :: ng
119    ! Cosine of solar zenith angle, single scattering albedo and
120    ! asymmetry factor:
121    real(jprb), intent(in)                :: mu0
122    real(jprb), intent(in),  dimension(:) :: ssa, g
123    real(jprb), intent(out), dimension(:) :: gamma1, gamma2, gamma3
124
125    real(jprb) :: factor
126
127    integer    :: jg
128
129#ifdef DO_DR_HOOK_TWO_STREAM
130    real(jphook) :: hook_handle
131
132    if (lhook) call dr_hook('radiation_two_stream:calc_two_stream_gammas_sw',0,hook_handle)
133#endif
134
135    ! Zdunkowski "PIFM" (Zdunkowski et al., 1980; Contributions to
136    ! Atmospheric Physics 53, 147-66)
137! Added for DWD (2020)
138!NEC$ shortloop
139    do jg = 1, ng
140      !      gamma1(jg) = 2.0_jprb  - ssa(jg) * (1.25_jprb + 0.75_jprb*g(jg))
141      !      gamma2(jg) = 0.75_jprb *(ssa(jg) * (1.0_jprb - g(jg)))
142      !      gamma3(jg) = 0.5_jprb  - (0.75_jprb*mu0)*g(jg)
143      ! Optimized version:
144      factor = 0.75_jprb*g(jg)
145      gamma1(jg) = 2.0_jprb  - ssa(jg) * (1.25_jprb + factor)
146      gamma2(jg) = ssa(jg) * (0.75_jprb - factor)
147      gamma3(jg) = 0.5_jprb  - mu0*factor
148    end do
149
150#ifdef DO_DR_HOOK_TWO_STREAM
151    if (lhook) call dr_hook('radiation_two_stream:calc_two_stream_gammas_sw',1,hook_handle)
152#endif
153
154  end subroutine calc_two_stream_gammas_sw
155
156
157  !---------------------------------------------------------------------
158  ! Compute the longwave reflectance and transmittance to diffuse
159  ! radiation using the Meador & Weaver formulas, as well as the
160  ! upward flux at the top and the downward flux at the base of the
161  ! layer due to emission from within the layer assuming a linear
162  ! variation of Planck function within the layer.
163  subroutine calc_reflectance_transmittance_lw(ng, &
164       &    od, gamma1, gamma2, planck_top, planck_bot, &
165       &    reflectance, transmittance, source_up, source_dn)
166
167#ifdef DO_DR_HOOK_TWO_STREAM
168    use yomhook, only : lhook, dr_hook, jphook
169#endif
170
171    implicit none
172   
173    integer, intent(in) :: ng
174
175    ! Optical depth and single scattering albedo
176    real(jprb), intent(in), dimension(ng) :: od
177
178    ! The two transfer coefficients from the two-stream
179    ! differentiatial equations (computed by
180    ! calc_two_stream_gammas_lw)
181    real(jprb), intent(in), dimension(ng) :: gamma1, gamma2
182
183    ! The Planck terms (functions of temperature) at the top and
184    ! bottom of the layer
185    real(jprb), intent(in), dimension(ng) :: planck_top, planck_bot
186
187    ! The diffuse reflectance and transmittance, i.e. the fraction of
188    ! diffuse radiation incident on a layer from either top or bottom
189    ! that is reflected back or transmitted through
190    real(jprb), intent(out), dimension(ng) :: reflectance, transmittance
191
192    ! The upward emission at the top of the layer and the downward
193    ! emission at its base, due to emission from within the layer
194    real(jprb), intent(out), dimension(ng) :: source_up, source_dn
195
196    real(jprd) :: k_exponent, reftrans_factor
197    real(jprd) :: exponential  ! = exp(-k_exponent*od)
198    real(jprd) :: exponential2 ! = exp(-2*k_exponent*od)
199
200    real(jprd) :: coeff, coeff_up_top, coeff_up_bot, coeff_dn_top, coeff_dn_bot
201
202    integer :: jg
203
204#ifdef DO_DR_HOOK_TWO_STREAM
205    real(jphook) :: hook_handle
206
207    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_lw',0,hook_handle)
208#endif
209
210! Added for DWD (2020)
211!NEC$ shortloop
212    do jg = 1, ng
213      if (od(jg) > 1.0e-3_jprd) then
214        k_exponent = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
215             1.0e-12_jprd)) ! Eq 18 of Meador & Weaver (1980)
216        exponential = exp_fast(-k_exponent*od(jg))
217        exponential2 = exponential*exponential
218        reftrans_factor = 1.0 / (k_exponent + gamma1(jg) + (k_exponent - gamma1(jg))*exponential2)
219        ! Meador & Weaver (1980) Eq. 25
220        reflectance(jg) = gamma2(jg) * (1.0_jprd - exponential2) * reftrans_factor
221        ! Meador & Weaver (1980) Eq. 26
222        transmittance(jg) = 2.0_jprd * k_exponent * exponential * reftrans_factor
223     
224        ! Compute upward and downward emission assuming the Planck
225        ! function to vary linearly with optical depth within the layer
226        ! (e.g. Wiscombe , JQSRT 1976).
227
228        ! Stackhouse and Stephens (JAS 1991) Eqs 5 & 12
229        coeff = (planck_bot(jg)-planck_top(jg)) / (od(jg)*(gamma1(jg)+gamma2(jg)))
230        coeff_up_top  =  coeff + planck_top(jg)
231        coeff_up_bot  =  coeff + planck_bot(jg)
232        coeff_dn_top  = -coeff + planck_top(jg)
233        coeff_dn_bot  = -coeff + planck_bot(jg)
234        source_up(jg) =  coeff_up_top - reflectance(jg) * coeff_dn_top - transmittance(jg) * coeff_up_bot
235        source_dn(jg) =  coeff_dn_bot - reflectance(jg) * coeff_up_bot - transmittance(jg) * coeff_dn_top
236      else
237        k_exponent = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
238             1.0e-12_jprd)) ! Eq 18 of Meador & Weaver (1980)
239        reflectance(jg) = gamma2(jg) * od(jg)
240        transmittance(jg) = (1.0_jprb - k_exponent*od(jg)) / (1.0_jprb + od(jg)*(gamma1(jg)-k_exponent))
241        source_up(jg) = (1.0_jprb - reflectance(jg) - transmittance(jg)) &
242             &       * 0.5 * (planck_top(jg) + planck_bot(jg))
243        source_dn(jg) = source_up(jg)
244      end if
245    end do
246   
247#ifdef DO_DR_HOOK_TWO_STREAM
248    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_lw',1,hook_handle)
249#endif
250 
251  end subroutine calc_reflectance_transmittance_lw
252 
253
254  !---------------------------------------------------------------------
255  ! Compute the longwave reflectance and transmittance to diffuse
256  ! radiation using the Meador & Weaver formulas, as well as the
257  ! upward flux at the top and the downward flux at the base of the
258  ! layer due to emission from within the layer assuming a linear
259  ! variation of Planck function within the layer; this version
260  ! computes gamma1 and gamma2 within the same routine.
261  subroutine calc_ref_trans_lw(ng, &
262       &    od, ssa, asymmetry, planck_top, planck_bot, &
263       &    reflectance, transmittance, source_up, source_dn)
264
265#ifdef DO_DR_HOOK_TWO_STREAM
266    use yomhook, only : lhook, dr_hook, jphook
267#endif
268
269    integer, intent(in) :: ng
270
271    ! Optical depth and single scattering albedo
272    real(jprb), intent(in), dimension(ng) :: od
273
274    ! Single scattering albedo and asymmetry factor
275    real(jprb), intent(in), dimension(ng) :: ssa, asymmetry
276
277    ! The Planck terms (functions of temperature) at the top and
278    ! bottom of the layer
279    real(jprb), intent(in), dimension(ng) :: planck_top, planck_bot
280
281    ! The diffuse reflectance and transmittance, i.e. the fraction of
282    ! diffuse radiation incident on a layer from either top or bottom
283    ! that is reflected back or transmitted through
284    real(jprb), intent(out), dimension(ng) :: reflectance, transmittance
285
286    ! The upward emission at the top of the layer and the downward
287    ! emission at its base, due to emission from within the layer
288    real(jprb), intent(out), dimension(ng) :: source_up, source_dn
289
290    ! The two transfer coefficients from the two-stream
291    ! differentiatial equations
292    real(jprb) :: gamma1, gamma2
293
294    real(jprb) :: k_exponent, reftrans_factor, factor
295    real(jprb) :: exponential  ! = exp(-k_exponent*od)
296    real(jprb) :: exponential2 ! = exp(-2*k_exponent*od)
297
298    real(jprb) :: coeff, coeff_up_top, coeff_up_bot, coeff_dn_top, coeff_dn_bot
299
300    integer :: jg
301
302#ifdef DO_DR_HOOK_TWO_STREAM
303    real(jphook) :: hook_handle
304
305    if (lhook) call dr_hook('radiation_two_stream:calc_ref_trans_lw',0,hook_handle)
306#endif
307
308! Added for DWD (2020)
309!NEC$ shortloop
310    do jg = 1, ng
311      factor = (LwDiffusivityWP * 0.5_jprb) * ssa(jg)
312      gamma1 = LwDiffusivityWP - factor*(1.0_jprb + asymmetry(jg))
313      gamma2 = factor * (1.0_jprb - asymmetry(jg))
314      k_exponent = sqrt(max((gamma1 - gamma2) * (gamma1 + gamma2), &
315           1.0e-12_jprb)) ! Eq 18 of Meador & Weaver (1980)
316      if (od(jg) > 1.0e-3_jprb) then
317        exponential = exp_fast(-k_exponent*od(jg))
318        exponential2 = exponential*exponential
319        reftrans_factor = 1.0_jprb / (k_exponent + gamma1 + (k_exponent - gamma1)*exponential2)
320        ! Meador & Weaver (1980) Eq. 25
321        reflectance(jg) = gamma2 * (1.0_jprb - exponential2) * reftrans_factor
322        ! Meador & Weaver (1980) Eq. 26
323        transmittance(jg) = 2.0_jprb * k_exponent * exponential * reftrans_factor
324     
325        ! Compute upward and downward emission assuming the Planck
326        ! function to vary linearly with optical depth within the layer
327        ! (e.g. Wiscombe , JQSRT 1976).
328
329        ! Stackhouse and Stephens (JAS 1991) Eqs 5 & 12
330        coeff = (planck_bot(jg)-planck_top(jg)) / (od(jg)*(gamma1+gamma2))
331        coeff_up_top  =  coeff + planck_top(jg)
332        coeff_up_bot  =  coeff + planck_bot(jg)
333        coeff_dn_top  = -coeff + planck_top(jg)
334        coeff_dn_bot  = -coeff + planck_bot(jg)
335        source_up(jg) =  coeff_up_top - reflectance(jg) * coeff_dn_top - transmittance(jg) * coeff_up_bot
336        source_dn(jg) =  coeff_dn_bot - reflectance(jg) * coeff_up_bot - transmittance(jg) * coeff_dn_top
337      else
338        reflectance(jg) = gamma2 * od(jg)
339        transmittance(jg) = (1.0_jprb - k_exponent*od(jg)) / (1.0_jprb + od(jg)*(gamma1-k_exponent))
340        source_up(jg) = (1.0_jprb - reflectance(jg) - transmittance(jg)) &
341             &       * 0.5 * (planck_top(jg) + planck_bot(jg))
342        source_dn(jg) = source_up(jg)
343      end if
344    end do
345   
346#ifdef DO_DR_HOOK_TWO_STREAM
347    if (lhook) call dr_hook('radiation_two_stream:calc_ref_trans_lw',1,hook_handle)
348#endif
349 
350  end subroutine calc_ref_trans_lw
351 
352 
353  !---------------------------------------------------------------------
354  ! Compute the longwave transmittance to diffuse radiation in the
355  ! no-scattering case, as well as the upward flux at the top and the
356  ! downward flux at the base of the layer due to emission from within
357  ! the layer assuming a linear variation of Planck function within
358  ! the layer.
359  subroutine calc_no_scattering_transmittance_lw(ng, &
360       &    od, planck_top, planck_bot, transmittance, source_up, source_dn)
361
362#ifdef DO_DR_HOOK_TWO_STREAM
363    use yomhook, only : lhook, dr_hook, jphook
364#endif
365
366    integer, intent(in) :: ng
367
368    ! Optical depth and single scattering albedo
369    real(jprb), intent(in), dimension(ng) :: od
370
371    ! The Planck terms (functions of temperature) at the top and
372    ! bottom of the layer
373    real(jprb), intent(in), dimension(ng) :: planck_top, planck_bot
374
375    ! The diffuse transmittance, i.e. the fraction of diffuse
376    ! radiation incident on a layer from either top or bottom that is
377    ! reflected back or transmitted through
378    real(jprb), intent(out), dimension(ng) :: transmittance
379
380    ! The upward emission at the top of the layer and the downward
381    ! emission at its base, due to emission from within the layer
382    real(jprb), intent(out), dimension(ng) :: source_up, source_dn
383
384    real(jprb) :: coeff, coeff_up_top, coeff_up_bot, coeff_dn_top, coeff_dn_bot !, planck_mean
385
386    integer :: jg
387
388#ifdef DO_DR_HOOK_TWO_STREAM
389    real(jphook) :: hook_handle
390
391    if (lhook) call dr_hook('radiation_two_stream:calc_no_scattering_transmittance_lw',0,hook_handle)
392#endif
393
394    transmittance = exp_fast(-LwDiffusivityWP*od)
395
396! Added for DWD (2020)
397!NEC$ shortloop
398    do jg = 1, ng
399      ! Compute upward and downward emission assuming the Planck
400      ! function to vary linearly with optical depth within the layer
401      ! (e.g. Wiscombe , JQSRT 1976).
402      coeff = LwDiffusivityWP*od(jg)
403      if (od(jg) > 1.0e-3_jprb) then
404        ! Simplified from calc_reflectance_transmittance_lw above
405        coeff = (planck_bot(jg)-planck_top(jg)) / coeff
406        coeff_up_top  =  coeff + planck_top(jg)
407        coeff_up_bot  =  coeff + planck_bot(jg)
408        coeff_dn_top  = -coeff + planck_top(jg)
409        coeff_dn_bot  = -coeff + planck_bot(jg)
410        source_up(jg) =  coeff_up_top - transmittance(jg) * coeff_up_bot
411        source_dn(jg) =  coeff_dn_bot - transmittance(jg) * coeff_dn_top
412      else
413        ! Linear limit at low optical depth
414        source_up(jg) = coeff * 0.5_jprb * (planck_top(jg)+planck_bot(jg))
415        source_dn(jg) = source_up(jg)
416      end if
417    end do
418
419#ifdef DO_DR_HOOK_TWO_STREAM
420    if (lhook) call dr_hook('radiation_two_stream:calc_no_scattering_transmittance_lw',1,hook_handle)
421#endif
422
423  end subroutine calc_no_scattering_transmittance_lw
424   
425   
426  !---------------------------------------------------------------------
427  ! Compute the shortwave reflectance and transmittance to diffuse
428  ! radiation using the Meador & Weaver formulas, as well as the
429  ! "direct" reflection and transmission, which really means the rate
430  ! of transfer of direct solar radiation (into a plane perpendicular
431  ! to the direct beam) into diffuse upward and downward streams at
432  ! the top and bottom of the layer, respectively.  Finally,
433  ! trans_dir_dir is the transmittance of the atmosphere to direct
434  ! radiation with no scattering.
435  subroutine calc_reflectance_transmittance_sw(ng, mu0, od, ssa, &
436       &      gamma1, gamma2, gamma3, ref_diff, trans_diff, &
437       &      ref_dir, trans_dir_diff, trans_dir_dir)
438   
439#ifdef DO_DR_HOOK_TWO_STREAM
440    use yomhook, only : lhook, dr_hook, jphook
441#endif
442
443    integer, intent(in) :: ng
444
445    ! Cosine of solar zenith angle
446    real(jprb), intent(in) :: mu0
447
448    ! Optical depth and single scattering albedo
449    real(jprb), intent(in), dimension(ng) :: od, ssa
450
451    ! The three transfer coefficients from the two-stream
452    ! differentiatial equations (computed by calc_two_stream_gammas)
453    real(jprb), intent(in), dimension(ng) :: gamma1, gamma2, gamma3
454
455    ! The direct reflectance and transmittance, i.e. the fraction of
456    ! incoming direct solar radiation incident at the top of a layer
457    ! that is either reflected back (ref_dir) or scattered but
458    ! transmitted through the layer to the base (trans_dir_diff)
459    real(jprb), intent(out), dimension(ng) :: ref_dir, trans_dir_diff
460
461    ! The diffuse reflectance and transmittance, i.e. the fraction of
462    ! diffuse radiation incident on a layer from either top or bottom
463    ! that is reflected back or transmitted through
464    real(jprb), intent(out), dimension(ng) :: ref_diff, trans_diff
465
466    ! Transmittance of the direct been with no scattering
467    real(jprb), intent(out), dimension(ng) :: trans_dir_dir
468
469    real(jprd) :: gamma4, alpha1, alpha2, k_exponent, reftrans_factor
470    real(jprb) :: exponential0 ! = exp(-od/mu0)
471    real(jprd) :: exponential  ! = exp(-k_exponent*od)
472    real(jprd) :: exponential2 ! = exp(-2*k_exponent*od)
473    real(jprd) :: k_mu0, k_gamma3, k_gamma4
474    real(jprd) :: k_2_exponential, od_over_mu0
475    integer    :: jg
476
477    ! Local value of cosine of solar zenith angle, in case it needs to be
478    ! tweaked to avoid near division by zero. This is intentionally in working
479    ! precision (jprb) rather than fixing at double precision (jprd).
480    real(jprb) :: mu0_local
481
482#ifdef DO_DR_HOOK_TWO_STREAM
483    real(jphook) :: hook_handle
484
485    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_sw',0,hook_handle)
486#endif
487
488! Added for DWD (2020)
489!NEC$ shortloop
490    do jg = 1, ng
491
492        gamma4 = 1.0_jprd - gamma3(jg)
493        alpha1 = gamma1(jg)*gamma4     + gamma2(jg)*gamma3(jg) ! Eq. 16
494        alpha2 = gamma1(jg)*gamma3(jg) + gamma2(jg)*gamma4    ! Eq. 17
495
496        k_exponent = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
497             &       1.0e-12_jprd)) ! Eq 18
498
499        ! We had a rare crash where k*mu0 was within around 1e-13 of 1,
500        ! leading to ref_dir and trans_dir_diff being well outside the range
501        ! 0-1. The following approach is appropriate when k_exponent is double
502        ! precision and mu0_local is single precision, although work is needed
503        ! to make this entire routine secure in single precision.
504        mu0_local = mu0
505        if (abs(1.0_jprd - k_exponent*mu0) < 1000.0_jprd * epsilon(1.0_jprd)) then
506          mu0_local = mu0 * (1.0_jprb - 10.0_jprb*epsilon(1.0_jprb))
507        end if
508
509        od_over_mu0 = max(od(jg) / mu0_local, 0.0_jprd)
510
511        ! Note that if the minimum value is reduced (e.g. to 1.0e-24)
512        ! then noise starts to appear as a function of solar zenith
513        ! angle
514        k_mu0 = k_exponent*mu0_local
515        k_gamma3 = k_exponent*gamma3(jg)
516        k_gamma4 = k_exponent*gamma4
517        ! Check for mu0 <= 0!
518        exponential0 = exp_fast(-od_over_mu0)
519        trans_dir_dir(jg) = exponential0
520        exponential = exp_fast(-k_exponent*od(jg))
521       
522        exponential2 = exponential*exponential
523        k_2_exponential = 2.0_jprd * k_exponent * exponential
524       
525        reftrans_factor = 1.0_jprd / (k_exponent + gamma1(jg) + (k_exponent - gamma1(jg))*exponential2)
526       
527        ! Meador & Weaver (1980) Eq. 25
528        ref_diff(jg) = gamma2(jg) * (1.0_jprd - exponential2) * reftrans_factor
529       
530        ! Meador & Weaver (1980) Eq. 26
531        trans_diff(jg) = k_2_exponential * reftrans_factor
532       
533        ! Here we need mu0 even though it wasn't in Meador and Weaver
534        ! because we are assuming the incoming direct flux is defined
535        ! to be the flux into a plane perpendicular to the direction of
536        ! the sun, not into a horizontal plane
537        reftrans_factor = mu0_local * ssa(jg) * reftrans_factor / (1.0_jprd - k_mu0*k_mu0)
538       
539        ! Meador & Weaver (1980) Eq. 14, multiplying top & bottom by
540        ! exp(-k_exponent*od) in case of very high optical depths
541        ref_dir(jg) = reftrans_factor &
542             &  * ( (1.0_jprd - k_mu0) * (alpha2 + k_gamma3) &
543             &     -(1.0_jprd + k_mu0) * (alpha2 - k_gamma3)*exponential2 &
544             &     -k_2_exponential*(gamma3(jg) - alpha2*mu0_local)*exponential0)
545       
546        ! Meador & Weaver (1980) Eq. 15, multiplying top & bottom by
547        ! exp(-k_exponent*od), minus the 1*exp(-od/mu0) term representing direct
548        ! unscattered transmittance. 
549        trans_dir_diff(jg) = reftrans_factor * ( k_2_exponential*(gamma4 + alpha1*mu0_local) &
550            & - exponential0 &
551            & * ( (1.0_jprd + k_mu0) * (alpha1 + k_gamma4) &
552            &    -(1.0_jprd - k_mu0) * (alpha1 - k_gamma4) * exponential2) )
553
554        ! Final check that ref_dir + trans_dir_diff <= 1
555        ref_dir(jg) = max(0.0_jprb, min(ref_dir(jg), 1.0_jprb))
556        trans_dir_diff(jg) = max(0.0_jprb, min(trans_dir_diff(jg), 1.0_jprb-ref_dir(jg)))
557
558    end do
559   
560#ifdef DO_DR_HOOK_TWO_STREAM
561    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_sw',1,hook_handle)
562#endif
563 
564  end subroutine calc_reflectance_transmittance_sw
565
566
567  !---------------------------------------------------------------------
568  ! Compute the shortwave reflectance and transmittance to diffuse
569  ! radiation using the Meador & Weaver formulas, as well as the
570  ! "direct" reflection and transmission, which really means the rate
571  ! of transfer of direct solar radiation (into a plane perpendicular
572  ! to the direct beam) into diffuse upward and downward streams at
573  ! the top and bottom of the layer, respectively.  Finally,
574  ! trans_dir_dir is the transmittance of the atmosphere to direct
575  ! radiation with no scattering. This version incorporates the
576  ! calculation of the gamma terms.
577  subroutine calc_ref_trans_sw(ng, mu0, od, ssa, &
578       &      asymmetry, ref_diff, trans_diff, &
579       &      ref_dir, trans_dir_diff, trans_dir_dir)
580   
581#ifdef DO_DR_HOOK_TWO_STREAM
582    use yomhook, only : lhook, dr_hook, jphook
583#endif
584
585    implicit none
586   
587    integer, intent(in) :: ng
588
589    ! Cosine of solar zenith angle
590    real(jprb), intent(in) :: mu0
591
592    ! Optical depth and single scattering albedo
593    real(jprb), intent(in), dimension(ng) :: od, ssa, asymmetry
594
595    ! The direct reflectance and transmittance, i.e. the fraction of
596    ! incoming direct solar radiation incident at the top of a layer
597    ! that is either reflected back (ref_dir) or scattered but
598    ! transmitted through the layer to the base (trans_dir_diff)
599    real(jprb), intent(out), dimension(ng) :: ref_dir, trans_dir_diff
600
601    ! The diffuse reflectance and transmittance, i.e. the fraction of
602    ! diffuse radiation incident on a layer from either top or bottom
603    ! that is reflected back or transmitted through
604    real(jprb), intent(out), dimension(ng) :: ref_diff, trans_diff
605
606    ! Transmittance of the direct been with no scattering
607    real(jprb), intent(out), dimension(ng) :: trans_dir_dir
608
609    ! The three transfer coefficients from the two-stream
610    ! differentiatial equations
611    real(jprb), dimension(ng) :: gamma1, gamma2, gamma3, gamma4
612    real(jprb), dimension(ng) :: alpha1, alpha2, k_exponent
613    real(jprb), dimension(ng) :: exponential ! = exp(-k_exponent*od)
614   
615    real(jprb) :: reftrans_factor, factor
616    real(jprb) :: exponential2 ! = exp(-2*k_exponent*od)
617    real(jprb) :: k_mu0, k_gamma3, k_gamma4
618    real(jprb) :: k_2_exponential, one_minus_kmu0_sqr
619    integer    :: jg
620
621#ifdef DO_DR_HOOK_TWO_STREAM
622    real(jphook) :: hook_handle
623
624    if (lhook) call dr_hook('radiation_two_stream:calc_ref_trans_sw',0,hook_handle)
625#endif
626
627    ! GCC 9.3 strange error: intermediate values of ~ -8000 cause a
628    ! FPE when vectorizing exp(), but not in non-vectorized loop, nor
629    ! with larger negative values!
630    trans_dir_dir = max(-max(od * (1.0_jprb/mu0), 0.0_jprb),-1000.0_jprb)
631    trans_dir_dir = exp_fast(trans_dir_dir)
632
633! Added for DWD (2020)
634!NEC$ shortloop
635    do jg = 1, ng
636
637      ! Zdunkowski "PIFM" (Zdunkowski et al., 1980; Contributions to
638      ! Atmospheric Physics 53, 147-66)
639      factor = 0.75_jprb*asymmetry(jg)
640
641      gamma1(jg) = 2.0_jprb  - ssa(jg) * (1.25_jprb + factor)
642      gamma2(jg) = ssa(jg) * (0.75_jprb - factor)
643      gamma3(jg) = 0.5_jprb  - mu0*factor
644      gamma4(jg) = 1.0_jprb - gamma3(jg)
645
646      alpha1(jg) = gamma1(jg)*gamma4(jg) + gamma2(jg)*gamma3(jg) ! Eq. 16
647      alpha2(jg) = gamma1(jg)*gamma3(jg) + gamma2(jg)*gamma4(jg) ! Eq. 17
648      ! The following line crashes inexplicably with gfortran 8.5.0 in
649      ! single precision - try a later version. Note that the minimum
650      ! value is needed to produce correct results for single
651      ! scattering albedos very close to or equal to one.
652#ifdef PARKIND1_SINGLE
653      k_exponent(jg) = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
654           &       1.0e-6_jprb)) ! Eq 18
655#else
656      k_exponent(jg) = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
657           &       1.0e-12_jprb)) ! Eq 18
658#endif
659    end do
660
661    exponential = exp_fast(-k_exponent*od)
662
663!NEC$ shortloop
664    do jg = 1, ng
665      k_mu0 = k_exponent(jg)*mu0
666      one_minus_kmu0_sqr = 1.0_jprb - k_mu0*k_mu0
667      k_gamma3 = k_exponent(jg)*gamma3(jg)
668      k_gamma4 = k_exponent(jg)*gamma4(jg)
669      exponential2 = exponential(jg)*exponential(jg)
670      k_2_exponential = 2.0_jprb * k_exponent(jg) * exponential(jg)
671      reftrans_factor = 1.0_jprb / (k_exponent(jg) + gamma1(jg) + (k_exponent(jg) - gamma1(jg))*exponential2)
672       
673      ! Meador & Weaver (1980) Eq. 25
674      ref_diff(jg) = gamma2(jg) * (1.0_jprb - exponential2) * reftrans_factor
675      !ref_diff(jg)       = max(0.0_jprb, min(ref_diff(jg)), 1.0_jprb)
676
677      ! Meador & Weaver (1980) Eq. 26, with security (which is
678      ! sometimes needed, but apparently not on ref_diff)
679      trans_diff(jg) = max(0.0_jprb, min(k_2_exponential * reftrans_factor, 1.0_jprb-ref_diff(jg)))
680
681      ! Here we need mu0 even though it wasn't in Meador and Weaver
682      ! because we are assuming the incoming direct flux is defined to
683      ! be the flux into a plane perpendicular to the direction of the
684      ! sun, not into a horizontal plane
685      reftrans_factor = mu0 * ssa(jg) * reftrans_factor &
686            &  / merge(one_minus_kmu0_sqr, epsilon(1.0_jprb), abs(one_minus_kmu0_sqr) > epsilon(1.0_jprb))
687     
688      ! Meador & Weaver (1980) Eq. 14, multiplying top & bottom by
689      ! exp(-k_exponent*od) in case of very high optical depths
690      ref_dir(jg) = reftrans_factor &
691           &  * ( (1.0_jprb - k_mu0) * (alpha2(jg) + k_gamma3) &
692           &     -(1.0_jprb + k_mu0) * (alpha2(jg) - k_gamma3)*exponential2 &
693           &     -k_2_exponential*(gamma3(jg) - alpha2(jg)*mu0)*trans_dir_dir(jg) )
694       
695      ! Meador & Weaver (1980) Eq. 15, multiplying top & bottom by
696      ! exp(-k_exponent*od), minus the 1*exp(-od/mu0) term
697      ! representing direct unscattered transmittance.
698      trans_dir_diff(jg) = reftrans_factor * ( k_2_exponential*(gamma4(jg) + alpha1(jg)*mu0) &
699           & - trans_dir_dir(jg) &
700           & * ( (1.0_jprb + k_mu0) * (alpha1(jg) + k_gamma4) &
701           &    -(1.0_jprb - k_mu0) * (alpha1(jg) - k_gamma4) * exponential2) )
702
703      ! Final check that ref_dir + trans_dir_diff <= 1
704      ref_dir(jg)        = max(0.0_jprb, min(ref_dir(jg), mu0*(1.0_jprb-trans_dir_dir(jg))))
705      trans_dir_diff(jg) = max(0.0_jprb, min(trans_dir_diff(jg), mu0*(1.0_jprb-trans_dir_dir(jg))-ref_dir(jg)))
706    end do
707   
708#ifdef DO_DR_HOOK_TWO_STREAM
709    if (lhook) call dr_hook('radiation_two_stream:calc_ref_trans_sw',1,hook_handle)
710#endif
711 
712  end subroutine calc_ref_trans_sw
713
714 
715  !---------------------------------------------------------------------
716  ! Compute the fraction of shortwave transmitted diffuse radiation
717  ! that is scattered during its transmission, used to compute
718  ! entrapment in SPARTACUS
719  subroutine calc_frac_scattered_diffuse_sw(ng, od, &
720       &      gamma1, gamma2, frac_scat_diffuse)
721   
722#ifdef DO_DR_HOOK_TWO_STREAM
723    use yomhook, only : lhook, dr_hook, jphook
724#endif
725
726    integer, intent(in) :: ng
727
728    ! Optical depth
729    real(jprb), intent(in), dimension(ng) :: od
730
731    ! The first two transfer coefficients from the two-stream
732    ! differentiatial equations (computed by calc_two_stream_gammas)
733    real(jprb), intent(in), dimension(ng) :: gamma1, gamma2
734
735    ! The fraction of shortwave transmitted diffuse radiation that is
736    ! scattered during its transmission
737    real(jprb), intent(out), dimension(ng) :: frac_scat_diffuse
738
739    real(jprd) :: k_exponent, reftrans_factor
740    real(jprd) :: exponential  ! = exp(-k_exponent*od)
741    real(jprd) :: exponential2 ! = exp(-2*k_exponent*od)
742    real(jprd) :: k_2_exponential
743    integer    :: jg
744
745#ifdef DO_DR_HOOK_TWO_STREAM
746    real(jphook) :: hook_handle
747
748    if (lhook) call dr_hook('radiation_two_stream:calc_frac_scattered_diffuse_sw',0,hook_handle)
749#endif
750
751! Added for DWD (2020)
752!NEC$ shortloop
753    do jg = 1, ng
754      ! Note that if the minimum value is reduced (e.g. to 1.0e-24)
755      ! then noise starts to appear as a function of solar zenith
756      ! angle
757      k_exponent = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
758           &       1.0e-12_jprd)) ! Eq 18
759      exponential = exp_fast(-k_exponent*od(jg))
760      exponential2 = exponential*exponential
761      k_2_exponential = 2.0_jprd * k_exponent * exponential
762       
763      reftrans_factor = 1.0_jprd / (k_exponent + gamma1(jg) + (k_exponent - gamma1(jg))*exponential2)
764       
765      ! Meador & Weaver (1980) Eq. 26.
766      ! Until 1.1.8, used LwDiffusivity instead of 2.0, although the
767      ! effect is very small
768      !      frac_scat_diffuse(jg) = 1.0_jprb - min(1.0_jprb,exp_fast(-LwDiffusivity*od(jg)) &
769      !           &  / max(1.0e-8_jprb, k_2_exponential * reftrans_factor))
770      frac_scat_diffuse(jg) = 1.0_jprb &
771           &  - min(1.0_jprb,exp_fast(-2.0_jprb*od(jg)) &
772           &  / max(1.0e-8_jprb, k_2_exponential * reftrans_factor))
773    end do
774   
775#ifdef DO_DR_HOOK_TWO_STREAM
776    if (lhook) call dr_hook('radiation_two_stream:calc_frac_scattered_diffuse_sw',1,hook_handle)
777#endif
778 
779  end subroutine calc_frac_scattered_diffuse_sw
780
781end module radiation_two_stream
Note: See TracBrowser for help on using the repository browser.