source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_liquid_optics_socrates.F90 @ 5445

Last change on this file since 5445 was 4773, checked in by idelkadi, 12 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: 2.9 KB
Line 
1! radiation_liquid_optics_socrates.F90 - SOCRATES method for parameterizing liquid droplet optics
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!   2020-08-10  R. Hogan  Bounded re to be >=1.2um and <=50um
17
18module radiation_liquid_optics_socrates
19
20  use parkind1, only : jprb
21
22  implicit none
23  public
24
25  ! SOCRATES (Edwards-Slingo) parameterizes info on the dependence of
26  ! the scattering properties in each band on effective radius in
27  ! terms of 16 coefficients
28  integer, parameter :: NLiqOpticsCoeffsSOCRATES = 16
29
30  ! Range of valid input effective radius, in microns
31  real(jprb), parameter :: MinEffectiveRadius = 1.2e-6
32  real(jprb), parameter :: MaxEffectiveRadius = 50.0e-6
33
34contains
35
36  !---------------------------------------------------------------------
37  ! Compute liquid-droplet scattering properties using a
38  ! parameterization consisting of Pade approximants from the
39  ! SOCRATES (Edwards-Slingo) code
40  subroutine calc_liq_optics_socrates(nb, coeff, lwp, re_in, od, scat_od, g)
41
42    use parkind1, only : jprb
43    !use yomhook,  only : lhook, dr_hook, jphook
44
45    ! Number of bands
46    integer, intent(in)  :: nb
47    ! Coefficients read from a data file
48    real(jprb), intent(in) :: coeff(:,:)
49    ! Liquid water path (kg m-2) and effective radius (m)
50    real(jprb), intent(in) :: lwp, re_in
51    ! Total optical depth, scattering optical depth and asymmetry factor
52    real(jprb), intent(out) :: od(nb), scat_od(nb), g(nb)
53
54    integer    :: jb
55    ! Local effective radius (m), after applying bounds
56    real(jprb) :: re
57
58    !real(jphook) :: hook_handle
59
60    !if (lhook) call dr_hook('radiation_liquid_optics_socrates:calc_liq_optics_socrates',0,hook_handle)
61
62    ! Apply the bounds of validity to effective radius
63    re = max(MinEffectiveRadius, min(re_in, MaxEffectiveRadius))
64
65! Added for DWD (2020)
66!NEC$ shortloop
67    do jb = 1, nb
68      od(jb) = lwp * (coeff(jb,1) + re*(coeff(jb,2) + re*coeff(jb,3))) &
69         &  / (1.0_jprb + re*(coeff(jb,4) + re*(coeff(jb,5) &
70         &  + re*coeff(jb,6))))
71      scat_od(jb) = od(jb) * (1.0_jprb &
72         &  - (coeff(jb,7) + re*(coeff(jb,8) + re*coeff(jb,9))) &
73         &  / (1.0_jprb + re*(coeff(jb,10) + re*coeff(jb,11))))
74      g(jb) = (coeff(jb,12) + re*(coeff(jb,13) + re*coeff(jb,14))) &
75         &  / (1.0_jprb + re*(coeff(jb,15) + re*coeff(jb,16)))
76    end do
77
78    !if (lhook) call dr_hook('radiation_liquid_optics_socrates:calc_liq_optics_socrates',1,hook_handle)
79
80  end subroutine calc_liq_optics_socrates
81
82end module radiation_liquid_optics_socrates
Note: See TracBrowser for help on using the repository browser.