source: LMDZ6/trunk/libf/phylmd/ecrad.v1.5.1/radiation_liquid_optics_slingo.F90 @ 4976

Last change on this file since 4976 was 3908, checked in by idelkadi, 3 years ago

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

  • Inclusion of the ecrad directory containing the sources of the ECRAD code
    • interface routine : radiation_scheme.F90
  • Adaptation of compilation scripts :
    • compilation under CPP key CPP_ECRAD
    • compilation with option "-rad ecard" or "-ecard true"
    • The "-rad old/rtm/ecran" build option will need to replace the "-rrtm true" and "-ecrad true" options in the future.
  • Runing LMDZ simulations with ecrad, you need :
    • logical key iflag_rrtm = 2 in physiq.def
    • namelist_ecrad (DefLists?)
    • the directory "data" containing the configuration files is temporarily placed in ../libfphylmd/ecrad/
  • Compilation and execution are tested in the 1D case. The repository under svn would allow to continue the implementation work: tests, verification of the results, ...
File size: 4.0 KB
Line 
1! radiation_liquid_optics_slingo.F90 - Slingo SW & Lindner-Li LW parameterization of liquid droplet optics
2!
3! (C) Copyright 2016- ECMWF.
4!
5! This software is licensed under the terms of the Apache Licence Version 2.0
6! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
7!
8! In applying this licence, ECMWF does not waive the privileges and immunities
9! granted to it by virtue of its status as an intergovernmental organisation
10! nor does it submit to any jurisdiction.
11!
12! Author:  Robin Hogan
13! Email:   r.j.hogan@ecmwf.int
14!
15
16module radiation_liquid_optics_slingo
17
18  implicit none
19  public
20
21  integer, parameter :: NLiqOpticsCoeffsSlingoSW = 6
22  integer, parameter :: NLiqOpticsCoeffsLindnerLiLW = 13
23
24contains
25
26  !---------------------------------------------------------------------
27  ! Compute liquid-droplet scattering properties in the shortwave from
28  ! Slingo (1989). WARNING: this parameterization is known not to be
29  ! very accurate: see Nielsen et al. (GMD 2014).
30  subroutine calc_liq_optics_slingo(nb, coeff, lwp, re, od, scat_od, g)
31
32    use parkind1, only : jprb
33    !use yomhook,  only : lhook, dr_hook
34
35    ! Number of bands
36    integer, intent(in)  :: nb
37    ! Coefficients read from a data file
38    real(jprb), intent(in) :: coeff(:,:)
39    ! Liquid water path (kg m-2) and effective radius (m)
40    real(jprb), intent(in) :: lwp, re
41    ! Total optical depth, scattering optical depth and asymmetry factor
42    real(jprb), intent(out) :: od(nb), scat_od(nb), g(nb)
43
44    ! Liquid water path in g m-2
45    real(jprb) :: lwp_gm_2
46    ! Effective radius in microns, and its inverse
47    real(jprb) :: re_um, inv_re_um
48
49    !real(jprb) :: hook_handle
50
51    !if (lhook) call dr_hook('radiation_liquid_optics_slingo:calc_liq_optics_slingo',0,hook_handle)
52
53    lwp_gm_2 = lwp * 1000.0_jprb
54    ! Range of validity reported by Slingo (1989): 4.2-16.6 microns
55    re_um = min(max(4.2_jprb, re * 1.0e6_jprb), 16.6_jprb)
56    inv_re_um = 1.0_jprb / re_um
57
58    od = lwp_gm_2 * (coeff(1:nb,1) + inv_re_um*coeff(1:nb,2))
59    scat_od = od * (1.0_jprb - coeff(1:nb,3) - re_um*coeff(1:nb,4))
60    g = coeff(1:nb,5) + re_um*coeff(1:nb,6)
61
62    !if (lhook) call dr_hook('radiation_liquid_optics_slingo:calc_liq_optics_slingo',1,hook_handle)
63
64  end subroutine calc_liq_optics_slingo
65
66
67  !---------------------------------------------------------------------
68  ! Compute liquid-droplet scattering properties in the longwave from
69  ! Lindner & Li (2000)
70  subroutine calc_liq_optics_lindner_li(nb, coeff, lwp, re, od, scat_od, g)
71
72    use parkind1, only : jprb
73    use yomhook,  only : lhook, dr_hook
74
75    ! Number of bands
76    integer, intent(in)  :: nb
77    ! Coefficients read from a data file
78    real(jprb), intent(in) :: coeff(:,:)
79    ! Liquid water path (kg m-2) and effective radius (m)
80    real(jprb), intent(in) :: lwp, re
81    ! Total optical depth, scattering optical depth and asymmetry factor
82    real(jprb), intent(out) :: od(nb), scat_od(nb), g(nb)
83
84    ! Liquid water path in g m-2
85    real(jprb) :: lwp_gm_2
86    ! Effective radius in microns, and its inverse
87    real(jprb) :: re_um, inv_re_um
88
89    real(jprb) :: hook_handle
90
91    if (lhook) call dr_hook('radiation_liquid_optics_slingo:calc_liq_optics_lindner_li',0,hook_handle)
92
93    lwp_gm_2 = lwp * 1000.0_jprb
94    ! Range of validity reported by Lindner and Li (2000): 2-40 microns
95    re_um = min(max(2.0_jprb, re * 1.0e6_jprb), 40.0_jprb)
96    inv_re_um = 1.0_jprb / re_um
97
98    od = lwp_gm_2 * (coeff(1:nb,1) + re_um*coeff(1:nb,2) + inv_re_um*(coeff(1:nb,3) &
99         &  + inv_re_um*(coeff(1:nb,4) + inv_re_um*coeff(1:nb,5))))
100    scat_od = od * (1.0_jprb - (coeff(1:nb,6) + inv_re_um*coeff(1:nb,7) &
101         &                      + re_um*(coeff(1:nb,8) + re_um*coeff(1:nb,9))))
102    g = coeff(1:nb,10) + inv_re_um*coeff(1:nb,11) &
103         &  + re_um*(coeff(1:nb,12) + re_um*coeff(1:nb,13))
104
105    if (lhook) call dr_hook('radiation_liquid_optics_slingo:calc_liq_optics_lindner_li',1,hook_handle)
106
107  end subroutine calc_liq_optics_lindner_li
108
109end module radiation_liquid_optics_slingo
Note: See TracBrowser for help on using the repository browser.