source: LMDZ6/trunk/libf/phylmd/ecrad.v1.5.1/radiation_two_stream.F90 @ 4947

Last change on this file since 4947 was 4489, checked in by lguez, 16 months ago

Merge LMDZ_ECRad branch back into trunk!

File size: 29.4 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
22module radiation_two_stream
23
24  use parkind1, only : jprb, jprd
25
26  implicit none
27  public
28
29  ! Elsasser's factor: the effective factor by which the zenith
30  ! optical depth needs to be multiplied to account for longwave
31  ! transmission at all angles through the atmosphere.  Alternatively
32  ! think of acos(1/lw_diffusivity) to be the effective zenith angle
33  ! of longwave radiation.
34  real(jprd), parameter :: LwDiffusivity   = 1.66_jprd
35  real(jprb), parameter :: LwDiffusivityWP = 1.66_jprb ! Working precision version
36
37  ! Shortwave diffusivity factor assumes hemispheric isotropy, assumed
38  ! by Zdunkowski's scheme and most others; note that for efficiency
39  ! this parameter is not used in the calculation of the gamma values,
40  ! but is used in the SPARTACUS solver.
41  real(jprb), parameter :: SwDiffusivity = 2.00_jprb
42
43  ! The routines in this module can be called millions of times, so
44  !calling Dr Hook for each one may be a significant overhead.
45  !Uncomment the following to turn Dr Hook on.
46!#define DO_DR_HOOK_TWO_STREAM
47
48contains
49
50#ifdef FAST_EXPONENTIAL
51  !---------------------------------------------------------------------
52  ! Fast exponential for negative arguments: a Pade approximant that
53  ! doesn't go negative for negative arguments, applied to arg/8, and
54  ! the result is then squared three times
55  elemental function exp_fast(arg) result(ex)
56    real(jprd) :: arg, ex
57    ex = 1.0_jprd / (1.0_jprd + arg*(-0.125_jprd &
58         + arg*(0.0078125_jprd - 0.000325520833333333_jprd * arg)))
59    ex = ex*ex
60    ex = ex*ex
61    ex = ex*ex
62  end function exp_fast
63#else
64#define exp_fast exp
65#endif
66
67  !---------------------------------------------------------------------
68  ! Calculate the two-stream coefficients gamma1 and gamma2 for the
69  ! longwave
70  subroutine calc_two_stream_gammas_lw(ng, ssa, g, &
71       &                               gamma1, gamma2)
72
73#ifdef DO_DR_HOOK_TWO_STREAM
74    use yomhook, only : lhook, dr_hook
75#endif
76
77    integer, intent(in) :: ng
78    ! Sngle scattering albedo and asymmetry factor:
79    real(jprb), intent(in),  dimension(:) :: ssa, g
80    real(jprb), intent(out), dimension(:) :: gamma1, gamma2
81
82    real(jprb) :: factor
83
84    integer    :: jg
85
86#ifdef DO_DR_HOOK_TWO_STREAM
87    real(jprb) :: hook_handle
88
89    if (lhook) call dr_hook('radiation_two_stream:calc_two_stream_gammas_lw',0,hook_handle)
90#endif
91! Added for DWD (2020)
92!NEC$ shortloop
93    do jg = 1, ng
94      ! Fu et al. (1997), Eq 2.9 and 2.10:
95      !      gamma1(jg) = LwDiffusivity * (1.0_jprb - 0.5_jprb*ssa(jg) &
96      !           &                    * (1.0_jprb + g(jg)))
97      !      gamma2(jg) = LwDiffusivity * 0.5_jprb * ssa(jg) &
98      !           &                    * (1.0_jprb - g(jg))
99      ! Reduce number of multiplications
100      factor = (LwDiffusivity * 0.5_jprb) * ssa(jg)
101      gamma1(jg) = LwDiffusivity - factor*(1.0_jprb + g(jg))
102      gamma2(jg) = factor * (1.0_jprb - g(jg))
103    end do
104
105#ifdef DO_DR_HOOK_TWO_STREAM
106    if (lhook) call dr_hook('radiation_two_stream:calc_two_stream_gammas_lw',1,hook_handle)
107#endif
108
109  end subroutine calc_two_stream_gammas_lw
110
111
112  !---------------------------------------------------------------------
113  ! Calculate the two-stream coefficients gamma1-gamma4 in the
114  ! shortwave
115  subroutine calc_two_stream_gammas_sw(ng, mu0, ssa, g, &
116       &                               gamma1, gamma2, gamma3)
117
118#ifdef DO_DR_HOOK_TWO_STREAM
119    use yomhook, only : lhook, dr_hook
120#endif
121
122    integer, intent(in) :: ng
123    ! Cosine of solar zenith angle, single scattering albedo and
124    ! asymmetry factor:
125    real(jprb), intent(in)                :: mu0
126    real(jprb), intent(in),  dimension(:) :: ssa, g
127    real(jprb), intent(out), dimension(:) :: gamma1, gamma2, gamma3
128
129    real(jprb) :: factor
130
131    integer    :: jg
132
133#ifdef DO_DR_HOOK_TWO_STREAM
134    real(jprb) :: hook_handle
135
136    if (lhook) call dr_hook('radiation_two_stream:calc_two_stream_gammas_sw',0,hook_handle)
137#endif
138
139    ! Zdunkowski "PIFM" (Zdunkowski et al., 1980; Contributions to
140    ! Atmospheric Physics 53, 147-66)
141! Added for DWD (2020)
142!NEC$ shortloop
143    do jg = 1, ng
144      !      gamma1(jg) = 2.0_jprb  - ssa(jg) * (1.25_jprb + 0.75_jprb*g(jg))
145      !      gamma2(jg) = 0.75_jprb *(ssa(jg) * (1.0_jprb - g(jg)))
146      !      gamma3(jg) = 0.5_jprb  - (0.75_jprb*mu0)*g(jg)
147      ! Optimized version:
148      factor = 0.75_jprb*g(jg)
149      gamma1(jg) = 2.0_jprb  - ssa(jg) * (1.25_jprb + factor)
150      gamma2(jg) = ssa(jg) * (0.75_jprb - factor)
151      gamma3(jg) = 0.5_jprb  - mu0*factor
152    end do
153
154#ifdef DO_DR_HOOK_TWO_STREAM
155    if (lhook) call dr_hook('radiation_two_stream:calc_two_stream_gammas_sw',1,hook_handle)
156#endif
157
158  end subroutine calc_two_stream_gammas_sw
159
160
161  !---------------------------------------------------------------------
162  ! Compute the longwave reflectance and transmittance to diffuse
163  ! radiation using the Meador & Weaver formulas, as well as the
164  ! upward flux at the top and the downward flux at the base of the
165  ! layer due to emission from within the layer assuming a linear
166  ! variation of Planck function within the layer.
167  subroutine calc_reflectance_transmittance_lw(ng, &
168       &    od, gamma1, gamma2, planck_top, planck_bot, &
169       &    reflectance, transmittance, source_up, source_dn)
170
171#ifdef DO_DR_HOOK_TWO_STREAM
172    use yomhook, only : lhook, dr_hook
173#endif
174
175    integer, intent(in) :: ng
176
177    ! Optical depth and single scattering albedo
178    real(jprb), intent(in), dimension(ng) :: od
179
180    ! The two transfer coefficients from the two-stream
181    ! differentiatial equations (computed by
182    ! calc_two_stream_gammas_lw)
183    real(jprb), intent(in), dimension(ng) :: gamma1, gamma2
184
185    ! The Planck terms (functions of temperature) at the top and
186    ! bottom of the layer
187    real(jprb), intent(in), dimension(ng) :: planck_top, planck_bot
188
189    ! The diffuse reflectance and transmittance, i.e. the fraction of
190    ! diffuse radiation incident on a layer from either top or bottom
191    ! that is reflected back or transmitted through
192    real(jprb), intent(out), dimension(ng) :: reflectance, transmittance
193
194    ! The upward emission at the top of the layer and the downward
195    ! emission at its base, due to emission from within the layer
196    real(jprb), intent(out), dimension(ng) :: source_up, source_dn
197
198    real(jprd) :: k_exponent, reftrans_factor
199    real(jprd) :: exponential  ! = exp(-k_exponent*od)
200    real(jprd) :: exponential2 ! = exp(-2*k_exponent*od)
201
202    real(jprd) :: coeff, coeff_up_top, coeff_up_bot, coeff_dn_top, coeff_dn_bot
203
204    integer :: jg
205
206#ifdef DO_DR_HOOK_TWO_STREAM
207    real(jprb) :: hook_handle
208
209    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_lw',0,hook_handle)
210#endif
211
212! Added for DWD (2020)
213!NEC$ shortloop
214    do jg = 1, ng
215      if (od(jg) > 1.0e-3_jprd) then
216        k_exponent = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
217             1.E-12_jprd)) ! Eq 18 of Meador & Weaver (1980)
218        exponential = exp_fast(-k_exponent*od(jg))
219        exponential2 = exponential*exponential
220        reftrans_factor = 1.0 / (k_exponent + gamma1(jg) + (k_exponent - gamma1(jg))*exponential2)
221        ! Meador & Weaver (1980) Eq. 25
222        reflectance(jg) = gamma2(jg) * (1.0_jprd - exponential2) * reftrans_factor
223        ! Meador & Weaver (1980) Eq. 26
224        transmittance(jg) = 2.0_jprd * k_exponent * exponential * reftrans_factor
225     
226        ! Compute upward and downward emission assuming the Planck
227        ! function to vary linearly with optical depth within the layer
228        ! (e.g. Wiscombe , JQSRT 1976).
229
230        ! Stackhouse and Stephens (JAS 1991) Eqs 5 & 12
231        coeff = (planck_bot(jg)-planck_top(jg)) / (od(jg)*(gamma1(jg)+gamma2(jg)))
232        coeff_up_top  =  coeff + planck_top(jg)
233        coeff_up_bot  =  coeff + planck_bot(jg)
234        coeff_dn_top  = -coeff + planck_top(jg)
235        coeff_dn_bot  = -coeff + planck_bot(jg)
236        source_up(jg) =  coeff_up_top - reflectance(jg) * coeff_dn_top - transmittance(jg) * coeff_up_bot
237        source_dn(jg) =  coeff_dn_bot - reflectance(jg) * coeff_up_bot - transmittance(jg) * coeff_dn_top
238      else
239        k_exponent = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
240             1.E-12_jprd)) ! Eq 18 of Meador & Weaver (1980)
241        reflectance(jg) = gamma2(jg) * od(jg)
242        transmittance(jg) = (1.0_jprb - k_exponent*od(jg)) / (1.0_jprb + od(jg)*(gamma1(jg)-k_exponent))
243        source_up(jg) = (1.0_jprb - reflectance(jg) - transmittance(jg)) &
244             &       * 0.5 * (planck_top(jg) + planck_bot(jg))
245        source_dn(jg) = source_up(jg)
246      end if
247    end do
248   
249#ifdef DO_DR_HOOK_TWO_STREAM
250    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_lw',1,hook_handle)
251#endif
252 
253  end subroutine calc_reflectance_transmittance_lw
254 
255
256
257  !---------------------------------------------------------------------
258  ! As calc_reflectance_transmittance_lw but for an isothermal layer
259  subroutine calc_reflectance_transmittance_isothermal_lw(ng, &
260       &    od, gamma1, gamma2, planck, &
261       &    reflectance, transmittance, source)
262
263#ifdef DO_DR_HOOK_TWO_STREAM
264    use yomhook, only : lhook, dr_hook
265#endif
266
267    integer, intent(in) :: ng
268
269    ! Optical depth and single scattering albedo
270    real(jprb), intent(in), dimension(ng) :: od
271
272    ! The two transfer coefficients from the two-stream
273    ! differentiatial equations (computed by
274    ! calc_two_stream_gammas_lw)
275    real(jprb), intent(in), dimension(ng) :: gamma1, gamma2
276
277    ! The Planck terms (functions of temperature) constant through the
278    ! layer
279    real(jprb), intent(in), dimension(ng) :: planck
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
289
290    real(jprd) :: k_exponent, reftrans_factor
291    real(jprd) :: exponential  ! = exp(-k_exponent*od)
292    real(jprd) :: exponential2 ! = exp(-2*k_exponent*od)
293
294    integer :: jg
295
296#ifdef DO_DR_HOOK_TWO_STREAM
297    real(jprb) :: hook_handle
298
299    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_isothermal_lw',0,hook_handle)
300#endif
301
302! Added for DWD (2020)
303!NEC$ shortloop
304    do jg = 1, ng
305      k_exponent = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
306           1.E-12_jprd)) ! Eq 18 of Meador & Weaver (1980)
307      exponential = exp_fast(-k_exponent*od(jg))
308      exponential2 = exponential*exponential
309      reftrans_factor = 1.0 / (k_exponent + gamma1(jg) + (k_exponent - gamma1(jg))*exponential2)
310      ! Meador & Weaver (1980) Eq. 25
311      reflectance(jg) = gamma2(jg) * (1.0_jprd - exponential2) * reftrans_factor
312      ! Meador & Weaver (1980) Eq. 26
313      transmittance(jg) = 2.0_jprd * k_exponent * exponential * reftrans_factor
314     
315      ! Emissivity of layer is one minus reflectance minus
316      ! transmittance, multiply by Planck function to get emitted
317      ! ousrce
318      source(jg) = planck(jg) * (1.0_jprd - reflectance(jg) - transmittance(jg))
319    end do
320   
321#ifdef DO_DR_HOOK_TWO_STREAM
322    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_isothermal_lw',1,hook_handle)
323#endif
324 
325  end subroutine calc_reflectance_transmittance_isothermal_lw
326 
327
328
329  !---------------------------------------------------------------------
330  ! Compute the longwave transmittance to diffuse radiation in the
331  ! no-scattering case, as well as the upward flux at the top and the
332  ! downward flux at the base of the layer due to emission from within
333  ! the layer assuming a linear variation of Planck function within
334  ! the layer.
335  subroutine calc_no_scattering_transmittance_lw(ng, &
336       &    od, planck_top, planck_bot, transmittance, source_up, source_dn)
337
338#ifdef DO_DR_HOOK_TWO_STREAM
339    use yomhook, only : lhook, dr_hook
340#endif
341
342    integer, intent(in) :: ng
343
344    ! Optical depth and single scattering albedo
345    real(jprb), intent(in), dimension(ng) :: od
346
347    ! The Planck terms (functions of temperature) at the top and
348    ! bottom of the layer
349    real(jprb), intent(in), dimension(ng) :: planck_top, planck_bot
350
351    ! The diffuse transmittance, i.e. the fraction of diffuse
352    ! radiation incident on a layer from either top or bottom that is
353    ! reflected back or transmitted through
354    real(jprb), intent(out), dimension(ng) :: transmittance
355
356    ! The upward emission at the top of the layer and the downward
357    ! emission at its base, due to emission from within the layer
358    real(jprb), intent(out), dimension(ng) :: source_up, source_dn
359
360    real(jprd) :: coeff, coeff_up_top, coeff_up_bot, coeff_dn_top, coeff_dn_bot !, planck_mean
361
362    integer :: jg
363
364#ifdef DO_DR_HOOK_TWO_STREAM
365    real(jprb) :: hook_handle
366
367    if (lhook) call dr_hook('radiation_two_stream:calc_no_scattering_transmittance_lw',0,hook_handle)
368#endif
369
370! Added for DWD (2020)
371!NEC$ shortloop
372    do jg = 1, ng
373      ! Compute upward and downward emission assuming the Planck
374      ! function to vary linearly with optical depth within the layer
375      ! (e.g. Wiscombe , JQSRT 1976).
376      if (od(jg) > 1.0e-3) then
377        ! Simplified from calc_reflectance_transmittance_lw above
378        coeff = LwDiffusivity*od(jg)
379        transmittance(jg) = exp_fast(-coeff)
380        coeff = (planck_bot(jg)-planck_top(jg)) / coeff
381        coeff_up_top  =  coeff + planck_top(jg)
382        coeff_up_bot  =  coeff + planck_bot(jg)
383        coeff_dn_top  = -coeff + planck_top(jg)
384        coeff_dn_bot  = -coeff + planck_bot(jg)
385        source_up(jg) =  coeff_up_top - transmittance(jg) * coeff_up_bot
386        source_dn(jg) =  coeff_dn_bot - transmittance(jg) * coeff_dn_top
387      else
388        ! Linear limit at low optical depth
389        coeff = LwDiffusivity*od(jg)
390        transmittance(jg) = 1.0_jprb - coeff
391        source_up(jg) = coeff * 0.5_jprb * (planck_top(jg)+planck_bot(jg))
392        source_dn(jg) = source_up(jg)
393      end if
394    end do
395
396    ! Method in the older IFS radiation scheme
397    !    do j = 1, n
398    !      coeff = od(jg) / (3.59712_jprd + od(jg))
399    !      planck_mean = 0.5_jprd * (planck_top(jg) + planck_bot(jg))
400    !     
401    !      source_up(jg) = (1.0_jprd-transmittance(jg)) * (planck_mean + (planck_top(jg)    - planck_mean) * coeff)
402    !      source_dn(jg) = (1.0_jprd-transmittance(jg)) * (planck_mean + (planck_bot(jg) - planck_mean) * coeff)
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
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(jprd) :: 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(jprb) :: 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_fast(-od_over_mu0)
505        trans_dir_dir(jg) = exponential0
506        exponential = exp_fast(-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  ! As above but with height as a vertical coordinate rather than
554  ! optical depth
555  subroutine calc_reflectance_transmittance_z_sw(ng, mu0, depth, &
556       &      gamma0, gamma1, gamma2, gamma3, gamma4, &
557       &      ref_diff, trans_diff, ref_dir, trans_dir_diff, trans_dir_dir)
558   
559#ifdef DO_DR_HOOK_TWO_STREAM
560    use yomhook, only : lhook, dr_hook
561#endif
562
563    integer, intent(in) :: ng
564
565    ! Cosine of solar zenith angle
566    real(jprb), intent(in) :: mu0
567
568    ! Layer depth
569    real(jprb), intent(in) :: depth
570
571    ! The four transfer coefficients from the two-stream
572    ! differentiatial equations
573    real(jprb), intent(in), dimension(ng) :: gamma1, gamma2, gamma3, gamma4
574
575    ! An additional coefficient for direct unscattered flux "Fdir"
576    ! such that dFdir/dz = -gamma0*Fdir
577    real(jprb), intent(in), dimension(ng) :: gamma0
578
579    ! The direct reflectance and transmittance, i.e. the fraction of
580    ! incoming direct solar radiation incident at the top of a layer
581    ! that is either reflected back (ref_dir) or scattered but
582    ! transmitted through the layer to the base (trans_dir_diff)
583    real(jprb), intent(out), dimension(ng) :: ref_dir, trans_dir_diff
584
585    ! The diffuse reflectance and transmittance, i.e. the fraction of
586    ! diffuse radiation incident on a layer from either top or bottom
587    ! that is reflected back or transmitted through
588    real(jprb), intent(out), dimension(ng) :: ref_diff, trans_diff
589
590    ! Transmittance of the direct been with no scattering
591    real(jprb), intent(out), dimension(ng) :: trans_dir_dir
592
593    real(jprd) :: alpha1, alpha2, k_exponent, reftrans_factor
594    real(jprd) :: exponential0 ! = exp(-od/mu0)
595    real(jprd) :: exponential  ! = exp(-k_exponent*od)
596    real(jprd) :: exponential2 ! = exp(-2*k_exponent*od)
597    real(jprd) :: k_mu0, k_gamma3, k_gamma4
598    real(jprd) :: k_2_exponential, od_over_mu0
599    integer    :: jg
600
601#ifdef DO_DR_HOOK_TWO_STREAM
602    real(jprb) :: hook_handle
603
604    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_z_sw',0,hook_handle)
605#endif
606
607! Added for DWD (2020)
608!NEC$ shortloop
609    do jg = 1, ng
610      od_over_mu0 = max(gamma0(jg) * depth, 0.0_jprd)
611      ! In the IFS this appears to be faster without testing the value
612      ! of od_over_mu0:
613      if (.true.) then
614!      if (od_over_mu0 > 1.0e-6_jprd) then
615        alpha1 = gamma1(jg)*gamma4(jg) + gamma2(jg)*gamma3(jg) ! Eq. 16
616        alpha2 = gamma1(jg)*gamma3(jg) + gamma2(jg)*gamma4(jg) ! Eq. 17
617       
618        ! Note that if the minimum value is reduced (e.g. to 1.0e-24)
619        ! then noise starts to appear as a function of solar zenith
620        ! angle
621        k_exponent = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
622             &       1.0e-12_jprd)) ! Eq 18
623        k_mu0 = k_exponent*mu0
624        k_gamma3 = k_exponent*gamma3(jg)
625        k_gamma4 = k_exponent*gamma4(jg)
626        ! Check for mu0 <= 0!
627        exponential0 = exp_fast(-od_over_mu0)
628        trans_dir_dir(jg) = exponential0
629        exponential = exp_fast(-k_exponent*depth)
630       
631        exponential2 = exponential*exponential
632        k_2_exponential = 2.0_jprd * k_exponent * exponential
633       
634        if (k_mu0 == 1.0_jprd) then
635          k_mu0 = 1.0_jprd - 10.0_jprd*epsilon(1.0_jprd)
636        end if
637       
638        reftrans_factor = 1.0_jprd / (k_exponent + gamma1(jg) + (k_exponent - gamma1(jg))*exponential2)
639       
640        ! Meador & Weaver (1980) Eq. 25
641        ref_diff(jg) = gamma2(jg) * (1.0_jprd - exponential2) * reftrans_factor
642       
643        ! Meador & Weaver (1980) Eq. 26
644        trans_diff(jg) = k_2_exponential * reftrans_factor
645       
646        ! Here we need mu0 even though it wasn't in Meador and Weaver
647        ! because we are assuming the incoming direct flux is defined
648        ! to be the flux into a plane perpendicular to the direction of
649        ! the sun, not into a horizontal plane
650        reftrans_factor = mu0 * reftrans_factor / (1.0_jprd - k_mu0*k_mu0)
651       
652        ! Meador & Weaver (1980) Eq. 14, multiplying top & bottom by
653        ! exp(-k_exponent*od) in case of very high optical depths
654        ref_dir(jg) = reftrans_factor &
655             &  * ( (1.0_jprd - k_mu0) * (alpha2 + k_gamma3) &
656             &     -(1.0_jprd + k_mu0) * (alpha2 - k_gamma3)*exponential2 &
657             &     -k_2_exponential*(gamma3(jg) - alpha2*mu0)*exponential0)
658       
659        ! Meador & Weaver (1980) Eq. 15, multiplying top & bottom by
660        ! exp(-k_exponent*od), minus the 1*exp(-od/mu0) term representing direct
661        ! unscattered transmittance. 
662        trans_dir_diff(jg) = reftrans_factor * ( k_2_exponential*(gamma4(jg) + alpha1*mu0) &
663            & - exponential0 &
664            & * ( (1.0_jprd + k_mu0) * (alpha1 + k_gamma4) &
665            &    -(1.0_jprd - k_mu0) * (alpha1 - k_gamma4) * exponential2) )
666
667      else
668        ! Low optical-depth limit; see equations 19, 20 and 27 from
669        ! Meador & Weaver (1980)
670        trans_diff(jg)     = 1.0_jprb - gamma1(jg) * depth
671        ref_diff(jg)       = gamma2(jg) * depth
672        trans_dir_diff(jg) = (1.0_jprb - gamma3(jg)) * depth
673        ref_dir(jg)        = gamma3(jg) * depth
674        trans_dir_dir(jg)  = 1.0_jprd - od_over_mu0
675      end if
676    end do
677   
678#ifdef DO_DR_HOOK_TWO_STREAM
679    if (lhook) call dr_hook('radiation_two_stream:calc_reflectance_transmittance_z_sw',1,hook_handle)
680#endif
681 
682  end subroutine calc_reflectance_transmittance_z_sw
683 
684
685  !---------------------------------------------------------------------
686  ! Compute the fraction of shortwave transmitted diffuse radiation
687  ! that is scattered during its transmission, used to compute
688  ! entrapment in SPARTACUS
689  subroutine calc_frac_scattered_diffuse_sw(ng, od, &
690       &      gamma1, gamma2, frac_scat_diffuse)
691   
692#ifdef DO_DR_HOOK_TWO_STREAM
693    use yomhook, only : lhook, dr_hook
694#endif
695
696    integer, intent(in) :: ng
697
698    ! Optical depth
699    real(jprb), intent(in), dimension(ng) :: od
700
701    ! The first two transfer coefficients from the two-stream
702    ! differentiatial equations (computed by calc_two_stream_gammas)
703    real(jprb), intent(in), dimension(ng) :: gamma1, gamma2
704
705    ! The fraction of shortwave transmitted diffuse radiation that is
706    ! scattered during its transmission
707    real(jprb), intent(out), dimension(ng) :: frac_scat_diffuse
708
709    real(jprd) :: k_exponent, reftrans_factor
710    real(jprd) :: exponential  ! = exp(-k_exponent*od)
711    real(jprd) :: exponential2 ! = exp(-2*k_exponent*od)
712    real(jprd) :: k_2_exponential
713    integer    :: jg
714
715#ifdef DO_DR_HOOK_TWO_STREAM
716    real(jprb) :: hook_handle
717
718    if (lhook) call dr_hook('radiation_two_stream:calc_frac_scattered_diffuse_sw',0,hook_handle)
719#endif
720
721! Added for DWD (2020)
722!NEC$ shortloop
723    do jg = 1, ng
724      ! Note that if the minimum value is reduced (e.g. to 1.0e-24)
725      ! then noise starts to appear as a function of solar zenith
726      ! angle
727      k_exponent = sqrt(max((gamma1(jg) - gamma2(jg)) * (gamma1(jg) + gamma2(jg)), &
728           &       1.0e-12_jprd)) ! Eq 18
729      exponential = exp_fast(-k_exponent*od(jg))
730      exponential2 = exponential*exponential
731      k_2_exponential = 2.0_jprd * k_exponent * exponential
732       
733      reftrans_factor = 1.0_jprd / (k_exponent + gamma1(jg) + (k_exponent - gamma1(jg))*exponential2)
734       
735      ! Meador & Weaver (1980) Eq. 26.
736      ! Until 1.1.8, used LwDiffusivity instead of 2.0, although the
737      ! effect is very small
738      !      frac_scat_diffuse(jg) = 1.0_jprb - min(1.0_jprb,exp_fast(-LwDiffusivity*od(jg)) &
739      !           &  / max(1.0e-8_jprb, k_2_exponential * reftrans_factor))
740      frac_scat_diffuse(jg) = 1.0_jprb &
741           &  - min(1.0_jprb,exp_fast(-2.0_jprb*od(jg)) &
742           &  / max(1.0e-8_jprb, k_2_exponential * reftrans_factor))
743    end do
744   
745#ifdef DO_DR_HOOK_TWO_STREAM
746    if (lhook) call dr_hook('radiation_two_stream:calc_frac_scattered_diffuse_sw',1,hook_handle)
747#endif
748 
749  end subroutine calc_frac_scattered_diffuse_sw
750
751end module radiation_two_stream
Note: See TracBrowser for help on using the repository browser.