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

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